xref: /petsc/src/ksp/pc/impls/asm/asm.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
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*/
1426cc229bSBarry Smith #include "petsc/private/matimpl.h"
154b9ad928SBarry Smith 
169371c9d4SSatish Balay static PetscErrorCode PCView_ASM(PC pc, PetscViewer viewer) {
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;
269566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
279566063dSJacob Faibussowitsch   PetscCall(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";
3063a3b9bcSJacob Faibussowitsch     if (osm->overlap >= 0) PetscCall(PetscSNPrintf(overlaps, sizeof(overlaps), "amount of overlap = %" PetscInt_FMT, osm->overlap));
3163a3b9bcSJacob Faibussowitsch     if (osm->n > 0) PetscCall(PetscSNPrintf(blocks, sizeof(blocks), "total subdomain blocks = %" PetscInt_FMT, osm->n));
329566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  %s, %s\n", blocks, overlaps));
339566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  restriction/interpolation type - %s\n", PCASMTypes[osm->type]));
349566063dSJacob Faibussowitsch     if (osm->dm_subdomains) PetscCall(PetscViewerASCIIPrintf(viewer, "  Additive Schwarz: using DM to define subdomains\n"));
359566063dSJacob Faibussowitsch     if (osm->loctype != PC_COMPOSITE_ADDITIVE) PetscCall(PetscViewerASCIIPrintf(viewer, "  Additive Schwarz: local solve composition type - %s\n", PCCompositeTypes[osm->loctype]));
369566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank));
379566063dSJacob Faibussowitsch     PetscCall(PetscViewerGetFormat(viewer, &format));
38ed155784SPierre Jolivet     if (format != PETSC_VIEWER_ASCII_INFO_DETAIL) {
394d219a6aSLisandro Dalcin       if (osm->ksp) {
409566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  Local solver information for first block is in the following KSP and PC objects on rank 0:\n"));
419566063dSJacob Faibussowitsch         PetscCall(PCGetOptionsPrefix(pc, &prefix));
429566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  Use -%sksp_view ::ascii_info_detail to display information for all blocks\n", prefix ? prefix : ""));
439566063dSJacob Faibussowitsch         PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
44dd400576SPatrick Sanan         if (rank == 0) {
459566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPushTab(viewer));
469566063dSJacob Faibussowitsch           PetscCall(KSPView(osm->ksp[0], sviewer));
479566063dSJacob Faibussowitsch           PetscCall(PetscViewerASCIIPopTab(viewer));
484b9ad928SBarry Smith         }
499566063dSJacob Faibussowitsch         PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
504d219a6aSLisandro Dalcin       }
514b9ad928SBarry Smith     } else {
529566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushSynchronized(viewer));
5363a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "  [%d] number of local blocks = %" PetscInt_FMT "\n", (int)rank, osm->n_local_true));
549566063dSJacob Faibussowitsch       PetscCall(PetscViewerFlush(viewer));
559566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  Local solver information for each block is in the following KSP and PC objects:\n"));
569566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushTab(viewer));
579566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "- - - - - - - - - - - - - - - - - -\n"));
589566063dSJacob Faibussowitsch       PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
59dd2fa690SBarry Smith       for (i = 0; i < osm->n_local_true; i++) {
609566063dSJacob Faibussowitsch         PetscCall(ISGetLocalSize(osm->is[i], &bsz));
6163a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIISynchronizedPrintf(sviewer, "[%d] local block number %" PetscInt_FMT ", size = %" PetscInt_FMT "\n", (int)rank, i, bsz));
629566063dSJacob Faibussowitsch         PetscCall(KSPView(osm->ksp[i], sviewer));
639566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIISynchronizedPrintf(sviewer, "- - - - - - - - - - - - - - - - - -\n"));
644b9ad928SBarry Smith       }
659566063dSJacob Faibussowitsch       PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
669566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopTab(viewer));
679566063dSJacob Faibussowitsch       PetscCall(PetscViewerFlush(viewer));
689566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopSynchronized(viewer));
694b9ad928SBarry Smith     }
704b9ad928SBarry Smith   } else if (isstring) {
7163a3b9bcSJacob Faibussowitsch     PetscCall(PetscViewerStringSPrintf(viewer, " blocks=%" PetscInt_FMT ", overlap=%" PetscInt_FMT ", type=%s", osm->n, osm->overlap, PCASMTypes[osm->type]));
729566063dSJacob Faibussowitsch     PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
739566063dSJacob Faibussowitsch     if (osm->ksp) PetscCall(KSPView(osm->ksp[0], sviewer));
749566063dSJacob Faibussowitsch     PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
754b9ad928SBarry Smith   }
764b9ad928SBarry Smith   PetscFunctionReturn(0);
774b9ad928SBarry Smith }
784b9ad928SBarry Smith 
799371c9d4SSatish Balay static PetscErrorCode PCASMPrintSubdomains(PC pc) {
8092bb6962SLisandro Dalcin   PC_ASM         *osm = (PC_ASM *)pc->data;
8192bb6962SLisandro Dalcin   const char     *prefix;
8292bb6962SLisandro Dalcin   char            fname[PETSC_MAX_PATH_LEN + 1];
83643535aeSDmitry Karpeev   PetscViewer     viewer, sviewer;
84643535aeSDmitry Karpeev   char           *s;
8592bb6962SLisandro Dalcin   PetscInt        i, j, nidx;
8692bb6962SLisandro Dalcin   const PetscInt *idx;
87643535aeSDmitry Karpeev   PetscMPIInt     rank, size;
88846783a0SBarry Smith 
8992bb6962SLisandro Dalcin   PetscFunctionBegin;
909566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size));
919566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank));
929566063dSJacob Faibussowitsch   PetscCall(PCGetOptionsPrefix(pc, &prefix));
939566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(NULL, prefix, "-pc_asm_print_subdomains", fname, sizeof(fname), NULL));
949566063dSJacob Faibussowitsch   if (fname[0] == 0) PetscCall(PetscStrcpy(fname, "stdout"));
959566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIOpen(PetscObjectComm((PetscObject)pc), fname, &viewer));
96643535aeSDmitry Karpeev   for (i = 0; i < osm->n_local; i++) {
97643535aeSDmitry Karpeev     if (i < osm->n_local_true) {
989566063dSJacob Faibussowitsch       PetscCall(ISGetLocalSize(osm->is[i], &nidx));
999566063dSJacob Faibussowitsch       PetscCall(ISGetIndices(osm->is[i], &idx));
100643535aeSDmitry Karpeev       /* Print to a string viewer; no more than 15 characters per index plus 512 char for the header.*/
10136a9e3b9SBarry Smith #define len 16 * (nidx + 1) + 512
1029566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(len, &s));
1039566063dSJacob Faibussowitsch       PetscCall(PetscViewerStringOpen(PETSC_COMM_SELF, s, len, &sviewer));
10436a9e3b9SBarry Smith #undef len
10563a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerStringSPrintf(sviewer, "[%d:%d] Subdomain %" PetscInt_FMT " with overlap:\n", rank, size, i));
106*48a46eb9SPierre Jolivet       for (j = 0; j < nidx; j++) PetscCall(PetscViewerStringSPrintf(sviewer, "%" PetscInt_FMT " ", idx[j]));
1079566063dSJacob Faibussowitsch       PetscCall(ISRestoreIndices(osm->is[i], &idx));
1089566063dSJacob Faibussowitsch       PetscCall(PetscViewerStringSPrintf(sviewer, "\n"));
1099566063dSJacob Faibussowitsch       PetscCall(PetscViewerDestroy(&sviewer));
1109566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushSynchronized(viewer));
11163a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "%s", s));
1129566063dSJacob Faibussowitsch       PetscCall(PetscViewerFlush(viewer));
1139566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopSynchronized(viewer));
1149566063dSJacob Faibussowitsch       PetscCall(PetscFree(s));
1152b691e39SMatthew Knepley       if (osm->is_local) {
116643535aeSDmitry Karpeev         /* Print to a string viewer; no more than 15 characters per index plus 512 char for the header.*/
11736a9e3b9SBarry Smith #define len 16 * (nidx + 1) + 512
1189566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(len, &s));
1199566063dSJacob Faibussowitsch         PetscCall(PetscViewerStringOpen(PETSC_COMM_SELF, s, len, &sviewer));
12036a9e3b9SBarry Smith #undef len
12163a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerStringSPrintf(sviewer, "[%d:%d] Subdomain %" PetscInt_FMT " without overlap:\n", rank, size, i));
1229566063dSJacob Faibussowitsch         PetscCall(ISGetLocalSize(osm->is_local[i], &nidx));
1239566063dSJacob Faibussowitsch         PetscCall(ISGetIndices(osm->is_local[i], &idx));
124*48a46eb9SPierre Jolivet         for (j = 0; j < nidx; j++) PetscCall(PetscViewerStringSPrintf(sviewer, "%" PetscInt_FMT " ", idx[j]));
1259566063dSJacob Faibussowitsch         PetscCall(ISRestoreIndices(osm->is_local[i], &idx));
1269566063dSJacob Faibussowitsch         PetscCall(PetscViewerStringSPrintf(sviewer, "\n"));
1279566063dSJacob Faibussowitsch         PetscCall(PetscViewerDestroy(&sviewer));
1289566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPushSynchronized(viewer));
12963a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "%s", s));
1309566063dSJacob Faibussowitsch         PetscCall(PetscViewerFlush(viewer));
1319566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPopSynchronized(viewer));
1329566063dSJacob Faibussowitsch         PetscCall(PetscFree(s));
133643535aeSDmitry Karpeev       }
1342fa5cd67SKarl Rupp     } else {
135643535aeSDmitry Karpeev       /* Participate in collective viewer calls. */
1369566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushSynchronized(viewer));
1379566063dSJacob Faibussowitsch       PetscCall(PetscViewerFlush(viewer));
1389566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopSynchronized(viewer));
139643535aeSDmitry Karpeev       /* Assume either all ranks have is_local or none do. */
140643535aeSDmitry Karpeev       if (osm->is_local) {
1419566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPushSynchronized(viewer));
1429566063dSJacob Faibussowitsch         PetscCall(PetscViewerFlush(viewer));
1439566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPopSynchronized(viewer));
144643535aeSDmitry Karpeev       }
145fdc48646SDmitry Karpeev     }
14692bb6962SLisandro Dalcin   }
1479566063dSJacob Faibussowitsch   PetscCall(PetscViewerFlush(viewer));
1489566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&viewer));
14992bb6962SLisandro Dalcin   PetscFunctionReturn(0);
15092bb6962SLisandro Dalcin }
15192bb6962SLisandro Dalcin 
1529371c9d4SSatish Balay static PetscErrorCode PCSetUp_ASM(PC pc) {
1534b9ad928SBarry Smith   PC_ASM     *osm = (PC_ASM *)pc->data;
15457501b6eSBarry Smith   PetscBool   flg;
15587e86712SBarry Smith   PetscInt    i, m, m_local;
1564b9ad928SBarry Smith   MatReuse    scall = MAT_REUSE_MATRIX;
1574b9ad928SBarry Smith   IS          isl;
1584b9ad928SBarry Smith   KSP         ksp;
1594b9ad928SBarry Smith   PC          subpc;
1602dcb1b2aSMatthew Knepley   const char *prefix, *pprefix;
16123ce1328SBarry Smith   Vec         vec;
1620298fd71SBarry Smith   DM         *domain_dm = NULL;
1634b9ad928SBarry Smith 
1644b9ad928SBarry Smith   PetscFunctionBegin;
1654b9ad928SBarry Smith   if (!pc->setupcalled) {
166265a2120SBarry Smith     PetscInt m;
16792bb6962SLisandro Dalcin 
1682b837212SDmitry Karpeev     /* Note: if subdomains have been set either via PCASMSetTotalSubdomains() or via PCASMSetLocalSubdomains(), osm->n_local_true will not be PETSC_DECIDE */
1692b837212SDmitry Karpeev     if (osm->n_local_true == PETSC_DECIDE) {
17069ca1f37SDmitry Karpeev       /* no subdomains given */
17165db9045SDmitry Karpeev       /* try pc->dm first, if allowed */
172d709ea83SDmitry Karpeev       if (osm->dm_subdomains && pc->dm) {
173feb221f8SDmitry Karpeev         PetscInt num_domains, d;
174feb221f8SDmitry Karpeev         char   **domain_names;
1758d4ac253SDmitry Karpeev         IS      *inner_domain_is, *outer_domain_is;
1769566063dSJacob Faibussowitsch         PetscCall(DMCreateDomainDecomposition(pc->dm, &num_domains, &domain_names, &inner_domain_is, &outer_domain_is, &domain_dm));
177704f0395SPatrick Sanan         osm->overlap = -1; /* We do not want to increase the overlap of the IS.
178704f0395SPatrick Sanan                               A future improvement of this code might allow one to use
179704f0395SPatrick Sanan                               DM-defined subdomains and also increase the overlap,
180704f0395SPatrick Sanan                               but that is not currently supported */
1811baa6e33SBarry Smith         if (num_domains) PetscCall(PCASMSetLocalSubdomains(pc, num_domains, outer_domain_is, inner_domain_is));
182feb221f8SDmitry Karpeev         for (d = 0; d < num_domains; ++d) {
1839566063dSJacob Faibussowitsch           if (domain_names) PetscCall(PetscFree(domain_names[d]));
1849566063dSJacob Faibussowitsch           if (inner_domain_is) PetscCall(ISDestroy(&inner_domain_is[d]));
1859566063dSJacob Faibussowitsch           if (outer_domain_is) PetscCall(ISDestroy(&outer_domain_is[d]));
186feb221f8SDmitry Karpeev         }
1879566063dSJacob Faibussowitsch         PetscCall(PetscFree(domain_names));
1889566063dSJacob Faibussowitsch         PetscCall(PetscFree(inner_domain_is));
1899566063dSJacob Faibussowitsch         PetscCall(PetscFree(outer_domain_is));
190feb221f8SDmitry Karpeev       }
1912b837212SDmitry Karpeev       if (osm->n_local_true == PETSC_DECIDE) {
19269ca1f37SDmitry Karpeev         /* still no subdomains; use one subdomain per processor */
1932b837212SDmitry Karpeev         osm->n_local_true = 1;
194feb221f8SDmitry Karpeev       }
1952b837212SDmitry Karpeev     }
1962b837212SDmitry Karpeev     { /* determine the global and max number of subdomains */
1979371c9d4SSatish Balay       struct {
1989371c9d4SSatish Balay         PetscInt max, sum;
1999371c9d4SSatish Balay       } inwork, outwork;
200c10200c1SHong Zhang       PetscMPIInt size;
201c10200c1SHong Zhang 
2026ac3741eSJed Brown       inwork.max = osm->n_local_true;
2036ac3741eSJed Brown       inwork.sum = osm->n_local_true;
2041c2dc1cbSBarry Smith       PetscCall(MPIU_Allreduce(&inwork, &outwork, 1, MPIU_2INT, MPIU_MAXSUM_OP, PetscObjectComm((PetscObject)pc)));
2056ac3741eSJed Brown       osm->n_local = outwork.max;
2066ac3741eSJed Brown       osm->n       = outwork.sum;
207c10200c1SHong Zhang 
2089566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size));
209c10200c1SHong Zhang       if (outwork.max == 1 && outwork.sum == size) {
2107dae84e0SHong Zhang         /* osm->n_local_true = 1 on all processes, set this option may enable use of optimized MatCreateSubMatrices() implementation */
2119566063dSJacob Faibussowitsch         PetscCall(MatSetOption(pc->pmat, MAT_SUBMAT_SINGLEIS, PETSC_TRUE));
212c10200c1SHong Zhang       }
2134b9ad928SBarry Smith     }
21478904715SLisandro Dalcin     if (!osm->is) { /* create the index sets */
2159566063dSJacob Faibussowitsch       PetscCall(PCASMCreateSubdomains(pc->pmat, osm->n_local_true, &osm->is));
2164b9ad928SBarry Smith     }
217f5234e35SJed Brown     if (osm->n_local_true > 1 && !osm->is_local) {
2189566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(osm->n_local_true, &osm->is_local));
219f5234e35SJed Brown       for (i = 0; i < osm->n_local_true; i++) {
220f5234e35SJed Brown         if (osm->overlap > 0) { /* With positive overlap, osm->is[i] will be modified */
2219566063dSJacob Faibussowitsch           PetscCall(ISDuplicate(osm->is[i], &osm->is_local[i]));
2229566063dSJacob Faibussowitsch           PetscCall(ISCopy(osm->is[i], osm->is_local[i]));
223f5234e35SJed Brown         } else {
2249566063dSJacob Faibussowitsch           PetscCall(PetscObjectReference((PetscObject)osm->is[i]));
225f5234e35SJed Brown           osm->is_local[i] = osm->is[i];
226f5234e35SJed Brown         }
227f5234e35SJed Brown       }
228f5234e35SJed Brown     }
2299566063dSJacob Faibussowitsch     PetscCall(PCGetOptionsPrefix(pc, &prefix));
2303d03bb48SJed Brown     if (osm->overlap > 0) {
2314b9ad928SBarry Smith       /* Extend the "overlapping" regions by a number of steps */
2329566063dSJacob Faibussowitsch       PetscCall(MatIncreaseOverlap(pc->pmat, osm->n_local_true, osm->is, osm->overlap));
2333d03bb48SJed Brown     }
2346ed231c7SMatthew Knepley     if (osm->sort_indices) {
23592bb6962SLisandro Dalcin       for (i = 0; i < osm->n_local_true; i++) {
2369566063dSJacob Faibussowitsch         PetscCall(ISSort(osm->is[i]));
237*48a46eb9SPierre Jolivet         if (osm->is_local) PetscCall(ISSort(osm->is_local[i]));
2384b9ad928SBarry Smith       }
2396ed231c7SMatthew Knepley     }
24098d3e228SPierre Jolivet     flg = PETSC_FALSE;
24198d3e228SPierre Jolivet     PetscCall(PetscOptionsHasName(NULL, prefix, "-pc_asm_print_subdomains", &flg));
24298d3e228SPierre Jolivet     if (flg) PetscCall(PCASMPrintSubdomains(pc));
243f6991133SBarry Smith     if (!osm->ksp) {
24478904715SLisandro Dalcin       /* Create the local solvers */
2459566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(osm->n_local_true, &osm->ksp));
246*48a46eb9SPierre Jolivet       if (domain_dm) PetscCall(PetscInfo(pc, "Setting up ASM subproblems using the embedded DM\n"));
24792bb6962SLisandro Dalcin       for (i = 0; i < osm->n_local_true; i++) {
2489566063dSJacob Faibussowitsch         PetscCall(KSPCreate(PETSC_COMM_SELF, &ksp));
2499566063dSJacob Faibussowitsch         PetscCall(KSPSetErrorIfNotConverged(ksp, pc->erroriffailure));
2509566063dSJacob Faibussowitsch         PetscCall(PetscLogObjectParent((PetscObject)pc, (PetscObject)ksp));
2519566063dSJacob Faibussowitsch         PetscCall(PetscObjectIncrementTabLevel((PetscObject)ksp, (PetscObject)pc, 1));
2529566063dSJacob Faibussowitsch         PetscCall(KSPSetType(ksp, KSPPREONLY));
2539566063dSJacob Faibussowitsch         PetscCall(KSPGetPC(ksp, &subpc));
2549566063dSJacob Faibussowitsch         PetscCall(PCGetOptionsPrefix(pc, &prefix));
2559566063dSJacob Faibussowitsch         PetscCall(KSPSetOptionsPrefix(ksp, prefix));
2569566063dSJacob Faibussowitsch         PetscCall(KSPAppendOptionsPrefix(ksp, "sub_"));
257feb221f8SDmitry Karpeev         if (domain_dm) {
2589566063dSJacob Faibussowitsch           PetscCall(KSPSetDM(ksp, domain_dm[i]));
2599566063dSJacob Faibussowitsch           PetscCall(KSPSetDMActive(ksp, PETSC_FALSE));
2609566063dSJacob Faibussowitsch           PetscCall(DMDestroy(&domain_dm[i]));
261feb221f8SDmitry Karpeev         }
2624b9ad928SBarry Smith         osm->ksp[i] = ksp;
2634b9ad928SBarry Smith       }
2641baa6e33SBarry Smith       if (domain_dm) PetscCall(PetscFree(domain_dm));
265f6991133SBarry Smith     }
2661dd8081eSeaulisa 
2679566063dSJacob Faibussowitsch     PetscCall(ISConcatenate(PETSC_COMM_SELF, osm->n_local_true, osm->is, &osm->lis));
2689566063dSJacob Faibussowitsch     PetscCall(ISSortRemoveDups(osm->lis));
2699566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(osm->lis, &m));
2701dd8081eSeaulisa 
2714b9ad928SBarry Smith     scall = MAT_INITIAL_MATRIX;
2724b9ad928SBarry Smith   } else {
2734b9ad928SBarry Smith     /*
2744b9ad928SBarry Smith        Destroy the blocks from the previous iteration
2754b9ad928SBarry Smith     */
276e09e08ccSBarry Smith     if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
2779566063dSJacob Faibussowitsch       PetscCall(MatDestroyMatrices(osm->n_local_true, &osm->pmat));
2784b9ad928SBarry Smith       scall = MAT_INITIAL_MATRIX;
2794b9ad928SBarry Smith     }
2804b9ad928SBarry Smith   }
2814b9ad928SBarry Smith 
282b58cb649SBarry Smith   /* Destroy previous submatrices of a different type than pc->pmat since MAT_REUSE_MATRIX won't work in that case */
28382b5ce2aSStefano Zampini   if (scall == MAT_REUSE_MATRIX && osm->sub_mat_type) {
284*48a46eb9SPierre Jolivet     if (osm->n_local_true > 0) PetscCall(MatDestroySubMatrices(osm->n_local_true, &osm->pmat));
285b58cb649SBarry Smith     scall = MAT_INITIAL_MATRIX;
286b58cb649SBarry Smith   }
287b58cb649SBarry Smith 
28892bb6962SLisandro Dalcin   /*
28992bb6962SLisandro Dalcin      Extract out the submatrices
29092bb6962SLisandro Dalcin   */
2919566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices(pc->pmat, osm->n_local_true, osm->is, osm->is, scall, &osm->pmat));
29292bb6962SLisandro Dalcin   if (scall == MAT_INITIAL_MATRIX) {
2939566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetOptionsPrefix((PetscObject)pc->pmat, &pprefix));
29492bb6962SLisandro Dalcin     for (i = 0; i < osm->n_local_true; i++) {
2959566063dSJacob Faibussowitsch       PetscCall(PetscLogObjectParent((PetscObject)pc, (PetscObject)osm->pmat[i]));
2969566063dSJacob Faibussowitsch       PetscCall(PetscObjectSetOptionsPrefix((PetscObject)osm->pmat[i], pprefix));
29792bb6962SLisandro Dalcin     }
29892bb6962SLisandro Dalcin   }
29980ec0b7dSPatrick Sanan 
30080ec0b7dSPatrick Sanan   /* Convert the types of the submatrices (if needbe) */
30180ec0b7dSPatrick Sanan   if (osm->sub_mat_type) {
302*48a46eb9SPierre Jolivet     for (i = 0; i < osm->n_local_true; i++) PetscCall(MatConvert(osm->pmat[i], osm->sub_mat_type, MAT_INPLACE_MATRIX, &(osm->pmat[i])));
30380ec0b7dSPatrick Sanan   }
30480ec0b7dSPatrick Sanan 
30580ec0b7dSPatrick Sanan   if (!pc->setupcalled) {
30656ea09ceSStefano Zampini     VecType vtype;
30756ea09ceSStefano Zampini 
30880ec0b7dSPatrick Sanan     /* Create the local work vectors (from the local matrices) and scatter contexts */
3099566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(pc->pmat, &vec, NULL));
3105b3c0d42Seaulisa 
3117827d75bSBarry Smith     PetscCheck(!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()");
312*48a46eb9SPierre Jolivet     if (osm->is_local && osm->type == PC_ASM_RESTRICT && osm->loctype == PC_COMPOSITE_ADDITIVE) PetscCall(PetscMalloc1(osm->n_local_true, &osm->lprolongation));
3139566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(osm->n_local_true, &osm->lrestriction));
3149566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(osm->n_local_true, &osm->x));
3159566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(osm->n_local_true, &osm->y));
316347574c9Seaulisa 
3179566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(osm->lis, &m));
3189566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_SELF, m, 0, 1, &isl));
3199566063dSJacob Faibussowitsch     PetscCall(MatGetVecType(osm->pmat[0], &vtype));
3209566063dSJacob Faibussowitsch     PetscCall(VecCreate(PETSC_COMM_SELF, &osm->lx));
3219566063dSJacob Faibussowitsch     PetscCall(VecSetSizes(osm->lx, m, m));
3229566063dSJacob Faibussowitsch     PetscCall(VecSetType(osm->lx, vtype));
3239566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(osm->lx, &osm->ly));
3249566063dSJacob Faibussowitsch     PetscCall(VecScatterCreate(vec, osm->lis, osm->lx, isl, &osm->restriction));
3259566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&isl));
326347574c9Seaulisa 
32780ec0b7dSPatrick Sanan     for (i = 0; i < osm->n_local_true; ++i) {
3285b3c0d42Seaulisa       ISLocalToGlobalMapping ltog;
3295b3c0d42Seaulisa       IS                     isll;
3305b3c0d42Seaulisa       const PetscInt        *idx_is;
3315b3c0d42Seaulisa       PetscInt              *idx_lis, nout;
3325b3c0d42Seaulisa 
3339566063dSJacob Faibussowitsch       PetscCall(ISGetLocalSize(osm->is[i], &m));
3349566063dSJacob Faibussowitsch       PetscCall(MatCreateVecs(osm->pmat[i], &osm->x[i], NULL));
3359566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(osm->x[i], &osm->y[i]));
33655817e79SBarry Smith 
337b0de9f23SBarry Smith       /* generate a scatter from ly to y[i] picking all the overlapping is[i] entries */
3389566063dSJacob Faibussowitsch       PetscCall(ISLocalToGlobalMappingCreateIS(osm->lis, &ltog));
3399566063dSJacob Faibussowitsch       PetscCall(ISGetLocalSize(osm->is[i], &m));
3409566063dSJacob Faibussowitsch       PetscCall(ISGetIndices(osm->is[i], &idx_is));
3419566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(m, &idx_lis));
3429566063dSJacob Faibussowitsch       PetscCall(ISGlobalToLocalMappingApply(ltog, IS_GTOLM_DROP, m, idx_is, &nout, idx_lis));
3437827d75bSBarry Smith       PetscCheck(nout == m, PETSC_COMM_SELF, PETSC_ERR_PLIB, "is not a subset of lis");
3449566063dSJacob Faibussowitsch       PetscCall(ISRestoreIndices(osm->is[i], &idx_is));
3459566063dSJacob Faibussowitsch       PetscCall(ISCreateGeneral(PETSC_COMM_SELF, m, idx_lis, PETSC_OWN_POINTER, &isll));
3469566063dSJacob Faibussowitsch       PetscCall(ISLocalToGlobalMappingDestroy(&ltog));
3479566063dSJacob Faibussowitsch       PetscCall(ISCreateStride(PETSC_COMM_SELF, m, 0, 1, &isl));
3489566063dSJacob Faibussowitsch       PetscCall(VecScatterCreate(osm->ly, isll, osm->y[i], isl, &osm->lrestriction[i]));
3499566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&isll));
3509566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&isl));
351910cf402Sprj-       if (osm->lprolongation) { /* generate a scatter from y[i] to ly picking only the the non-overlapping is_local[i] entries */
352d8b3b5e3Seaulisa         ISLocalToGlobalMapping ltog;
353d8b3b5e3Seaulisa         IS                     isll, isll_local;
354d8b3b5e3Seaulisa         const PetscInt        *idx_local;
355d8b3b5e3Seaulisa         PetscInt              *idx1, *idx2, nout;
356d8b3b5e3Seaulisa 
3579566063dSJacob Faibussowitsch         PetscCall(ISGetLocalSize(osm->is_local[i], &m_local));
3589566063dSJacob Faibussowitsch         PetscCall(ISGetIndices(osm->is_local[i], &idx_local));
359d8b3b5e3Seaulisa 
3609566063dSJacob Faibussowitsch         PetscCall(ISLocalToGlobalMappingCreateIS(osm->is[i], &ltog));
3619566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(m_local, &idx1));
3629566063dSJacob Faibussowitsch         PetscCall(ISGlobalToLocalMappingApply(ltog, IS_GTOLM_DROP, m_local, idx_local, &nout, idx1));
3639566063dSJacob Faibussowitsch         PetscCall(ISLocalToGlobalMappingDestroy(&ltog));
3647827d75bSBarry Smith         PetscCheck(nout == m_local, PETSC_COMM_SELF, PETSC_ERR_PLIB, "is_local not a subset of is");
3659566063dSJacob Faibussowitsch         PetscCall(ISCreateGeneral(PETSC_COMM_SELF, m_local, idx1, PETSC_OWN_POINTER, &isll));
366d8b3b5e3Seaulisa 
3679566063dSJacob Faibussowitsch         PetscCall(ISLocalToGlobalMappingCreateIS(osm->lis, &ltog));
3689566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(m_local, &idx2));
3699566063dSJacob Faibussowitsch         PetscCall(ISGlobalToLocalMappingApply(ltog, IS_GTOLM_DROP, m_local, idx_local, &nout, idx2));
3709566063dSJacob Faibussowitsch         PetscCall(ISLocalToGlobalMappingDestroy(&ltog));
3717827d75bSBarry Smith         PetscCheck(nout == m_local, PETSC_COMM_SELF, PETSC_ERR_PLIB, "is_local not a subset of lis");
3729566063dSJacob Faibussowitsch         PetscCall(ISCreateGeneral(PETSC_COMM_SELF, m_local, idx2, PETSC_OWN_POINTER, &isll_local));
373d8b3b5e3Seaulisa 
3749566063dSJacob Faibussowitsch         PetscCall(ISRestoreIndices(osm->is_local[i], &idx_local));
3759566063dSJacob Faibussowitsch         PetscCall(VecScatterCreate(osm->y[i], isll, osm->ly, isll_local, &osm->lprolongation[i]));
376d8b3b5e3Seaulisa 
3779566063dSJacob Faibussowitsch         PetscCall(ISDestroy(&isll));
3789566063dSJacob Faibussowitsch         PetscCall(ISDestroy(&isll_local));
379d8b3b5e3Seaulisa       }
38080ec0b7dSPatrick Sanan     }
3819566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&vec));
38280ec0b7dSPatrick Sanan   }
38380ec0b7dSPatrick Sanan 
384fb745f2cSMatthew G. Knepley   if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE) {
385235cc4ceSMatthew G. Knepley     IS      *cis;
386235cc4ceSMatthew G. Knepley     PetscInt c;
387235cc4ceSMatthew G. Knepley 
3889566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(osm->n_local_true, &cis));
389235cc4ceSMatthew G. Knepley     for (c = 0; c < osm->n_local_true; ++c) cis[c] = osm->lis;
3909566063dSJacob Faibussowitsch     PetscCall(MatCreateSubMatrices(pc->pmat, osm->n_local_true, osm->is, cis, scall, &osm->lmats));
3919566063dSJacob Faibussowitsch     PetscCall(PetscFree(cis));
392fb745f2cSMatthew G. Knepley   }
3934b9ad928SBarry Smith 
3944b9ad928SBarry Smith   /* Return control to the user so that the submatrices can be modified (e.g., to apply
3954b9ad928SBarry Smith      different boundary conditions for the submatrices than for the global problem) */
3969566063dSJacob Faibussowitsch   PetscCall(PCModifySubMatrices(pc, osm->n_local_true, osm->is, osm->is, osm->pmat, pc->modifysubmatricesP));
3974b9ad928SBarry Smith 
39892bb6962SLisandro Dalcin   /*
39992bb6962SLisandro Dalcin      Loop over subdomains putting them into local ksp
40092bb6962SLisandro Dalcin   */
4019566063dSJacob Faibussowitsch   PetscCall(KSPGetOptionsPrefix(osm->ksp[0], &prefix));
40292bb6962SLisandro Dalcin   for (i = 0; i < osm->n_local_true; i++) {
4039566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(osm->ksp[i], osm->pmat[i], osm->pmat[i]));
4049566063dSJacob Faibussowitsch     PetscCall(MatSetOptionsPrefix(osm->pmat[i], prefix));
405*48a46eb9SPierre Jolivet     if (!pc->setupcalled) PetscCall(KSPSetFromOptions(osm->ksp[i]));
40692bb6962SLisandro Dalcin   }
4074b9ad928SBarry Smith   PetscFunctionReturn(0);
4084b9ad928SBarry Smith }
4094b9ad928SBarry Smith 
4109371c9d4SSatish Balay static PetscErrorCode PCSetUpOnBlocks_ASM(PC pc) {
4114b9ad928SBarry Smith   PC_ASM            *osm = (PC_ASM *)pc->data;
41213f74950SBarry Smith   PetscInt           i;
413539c167fSBarry Smith   KSPConvergedReason reason;
4144b9ad928SBarry Smith 
4154b9ad928SBarry Smith   PetscFunctionBegin;
4164b9ad928SBarry Smith   for (i = 0; i < osm->n_local_true; i++) {
4179566063dSJacob Faibussowitsch     PetscCall(KSPSetUp(osm->ksp[i]));
4189566063dSJacob Faibussowitsch     PetscCall(KSPGetConvergedReason(osm->ksp[i], &reason));
4199371c9d4SSatish Balay     if (reason == KSP_DIVERGED_PC_FAILED) { pc->failedreason = PC_SUBPC_ERROR; }
4204b9ad928SBarry Smith   }
4214b9ad928SBarry Smith   PetscFunctionReturn(0);
4224b9ad928SBarry Smith }
4234b9ad928SBarry Smith 
4249371c9d4SSatish Balay static PetscErrorCode PCApply_ASM(PC pc, Vec x, Vec y) {
4254b9ad928SBarry Smith   PC_ASM     *osm = (PC_ASM *)pc->data;
4261dd8081eSeaulisa   PetscInt    i, n_local_true = osm->n_local_true;
4274b9ad928SBarry Smith   ScatterMode forward = SCATTER_FORWARD, reverse = SCATTER_REVERSE;
4284b9ad928SBarry Smith 
4294b9ad928SBarry Smith   PetscFunctionBegin;
4304b9ad928SBarry Smith   /*
43148e38a8aSPierre Jolivet      support for limiting the restriction or interpolation to only local
4324b9ad928SBarry Smith      subdomain values (leaving the other values 0).
4334b9ad928SBarry Smith   */
4344b9ad928SBarry Smith   if (!(osm->type & PC_ASM_RESTRICT)) {
4354b9ad928SBarry Smith     forward = SCATTER_FORWARD_LOCAL;
4364b9ad928SBarry Smith     /* have to zero the work RHS since scatter may leave some slots empty */
4379566063dSJacob Faibussowitsch     PetscCall(VecSet(osm->lx, 0.0));
4384b9ad928SBarry Smith   }
4399371c9d4SSatish Balay   if (!(osm->type & PC_ASM_INTERPOLATE)) { reverse = SCATTER_REVERSE_LOCAL; }
4404b9ad928SBarry Smith 
441347574c9Seaulisa   if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE || osm->loctype == PC_COMPOSITE_ADDITIVE) {
442b0de9f23SBarry Smith     /* zero the global and the local solutions */
4439566063dSJacob Faibussowitsch     PetscCall(VecSet(y, 0.0));
4449566063dSJacob Faibussowitsch     PetscCall(VecSet(osm->ly, 0.0));
445347574c9Seaulisa 
44648e38a8aSPierre Jolivet     /* copy the global RHS to local RHS including the ghost nodes */
4479566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(osm->restriction, x, osm->lx, INSERT_VALUES, forward));
4489566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(osm->restriction, x, osm->lx, INSERT_VALUES, forward));
449347574c9Seaulisa 
45048e38a8aSPierre Jolivet     /* restrict local RHS to the overlapping 0-block RHS */
4519566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward));
4529566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward));
453d8b3b5e3Seaulisa 
45412cd4985SMatthew G. Knepley     /* do the local solves */
45512cd4985SMatthew G. Knepley     for (i = 0; i < n_local_true; ++i) {
456b0de9f23SBarry Smith       /* solve the overlapping i-block */
4579566063dSJacob Faibussowitsch       PetscCall(PetscLogEventBegin(PC_ApplyOnBlocks, osm->ksp[i], osm->x[i], osm->y[i], 0));
4589566063dSJacob Faibussowitsch       PetscCall(KSPSolve(osm->ksp[i], osm->x[i], osm->y[i]));
4599566063dSJacob Faibussowitsch       PetscCall(KSPCheckSolve(osm->ksp[i], pc, osm->y[i]));
4609566063dSJacob Faibussowitsch       PetscCall(PetscLogEventEnd(PC_ApplyOnBlocks, osm->ksp[i], osm->x[i], osm->y[i], 0));
461d8b3b5e3Seaulisa 
462910cf402Sprj-       if (osm->lprolongation) { /* interpolate the non-overlapping i-block solution to the local solution (only for restrictive additive) */
4639566063dSJacob Faibussowitsch         PetscCall(VecScatterBegin(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward));
4649566063dSJacob Faibussowitsch         PetscCall(VecScatterEnd(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward));
46548e38a8aSPierre Jolivet       } else { /* interpolate the overlapping i-block solution to the local solution */
4669566063dSJacob Faibussowitsch         PetscCall(VecScatterBegin(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse));
4679566063dSJacob Faibussowitsch         PetscCall(VecScatterEnd(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse));
468d8b3b5e3Seaulisa       }
469347574c9Seaulisa 
470347574c9Seaulisa       if (i < n_local_true - 1) {
47148e38a8aSPierre Jolivet         /* restrict local RHS to the overlapping (i+1)-block RHS */
4729566063dSJacob Faibussowitsch         PetscCall(VecScatterBegin(osm->lrestriction[i + 1], osm->lx, osm->x[i + 1], INSERT_VALUES, forward));
4739566063dSJacob Faibussowitsch         PetscCall(VecScatterEnd(osm->lrestriction[i + 1], osm->lx, osm->x[i + 1], INSERT_VALUES, forward));
474347574c9Seaulisa 
475347574c9Seaulisa         if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE) {
4768966356dSPierre Jolivet           /* update the overlapping (i+1)-block RHS using the current local solution */
4779566063dSJacob Faibussowitsch           PetscCall(MatMult(osm->lmats[i + 1], osm->ly, osm->y[i + 1]));
4789566063dSJacob Faibussowitsch           PetscCall(VecAXPBY(osm->x[i + 1], -1., 1., osm->y[i + 1]));
4797c3d802fSMatthew G. Knepley         }
48012cd4985SMatthew G. Knepley       }
48112cd4985SMatthew G. Knepley     }
48248e38a8aSPierre Jolivet     /* add the local solution to the global solution including the ghost nodes */
4839566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(osm->restriction, osm->ly, y, ADD_VALUES, reverse));
4849566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(osm->restriction, osm->ly, y, ADD_VALUES, reverse));
48598921bdaSJacob Faibussowitsch   } else SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Invalid local composition type: %s", PCCompositeTypes[osm->loctype]);
4864b9ad928SBarry Smith   PetscFunctionReturn(0);
4874b9ad928SBarry Smith }
4884b9ad928SBarry Smith 
4899371c9d4SSatish Balay static PetscErrorCode PCMatApply_ASM(PC pc, Mat X, Mat Y) {
49048e38a8aSPierre Jolivet   PC_ASM     *osm = (PC_ASM *)pc->data;
49148e38a8aSPierre Jolivet   Mat         Z, W;
49248e38a8aSPierre Jolivet   Vec         x;
49348e38a8aSPierre Jolivet   PetscInt    i, m, N;
49448e38a8aSPierre Jolivet   ScatterMode forward = SCATTER_FORWARD, reverse = SCATTER_REVERSE;
49548e38a8aSPierre Jolivet 
49648e38a8aSPierre Jolivet   PetscFunctionBegin;
4977827d75bSBarry Smith   PetscCheck(osm->n_local_true <= 1, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Not yet implemented");
49848e38a8aSPierre Jolivet   /*
49948e38a8aSPierre Jolivet      support for limiting the restriction or interpolation to only local
50048e38a8aSPierre Jolivet      subdomain values (leaving the other values 0).
50148e38a8aSPierre Jolivet   */
50248e38a8aSPierre Jolivet   if (!(osm->type & PC_ASM_RESTRICT)) {
50348e38a8aSPierre Jolivet     forward = SCATTER_FORWARD_LOCAL;
50448e38a8aSPierre Jolivet     /* have to zero the work RHS since scatter may leave some slots empty */
5059566063dSJacob Faibussowitsch     PetscCall(VecSet(osm->lx, 0.0));
50648e38a8aSPierre Jolivet   }
5079371c9d4SSatish Balay   if (!(osm->type & PC_ASM_INTERPOLATE)) { reverse = SCATTER_REVERSE_LOCAL; }
5089566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(osm->x[0], &m));
5099566063dSJacob Faibussowitsch   PetscCall(MatGetSize(X, NULL, &N));
5109566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, m, N, NULL, &Z));
51148e38a8aSPierre Jolivet   if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE || osm->loctype == PC_COMPOSITE_ADDITIVE) {
51248e38a8aSPierre Jolivet     /* zero the global and the local solutions */
5139566063dSJacob Faibussowitsch     PetscCall(MatZeroEntries(Y));
5149566063dSJacob Faibussowitsch     PetscCall(VecSet(osm->ly, 0.0));
51548e38a8aSPierre Jolivet 
51648e38a8aSPierre Jolivet     for (i = 0; i < N; ++i) {
5179566063dSJacob Faibussowitsch       PetscCall(MatDenseGetColumnVecRead(X, i, &x));
51848e38a8aSPierre Jolivet       /* copy the global RHS to local RHS including the ghost nodes */
5199566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(osm->restriction, x, osm->lx, INSERT_VALUES, forward));
5209566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(osm->restriction, x, osm->lx, INSERT_VALUES, forward));
5219566063dSJacob Faibussowitsch       PetscCall(MatDenseRestoreColumnVecRead(X, i, &x));
52248e38a8aSPierre Jolivet 
5239566063dSJacob Faibussowitsch       PetscCall(MatDenseGetColumnVecWrite(Z, i, &x));
52448e38a8aSPierre Jolivet       /* restrict local RHS to the overlapping 0-block RHS */
5259566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(osm->lrestriction[0], osm->lx, x, INSERT_VALUES, forward));
5269566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(osm->lrestriction[0], osm->lx, x, INSERT_VALUES, forward));
5279566063dSJacob Faibussowitsch       PetscCall(MatDenseRestoreColumnVecWrite(Z, i, &x));
52848e38a8aSPierre Jolivet     }
5299566063dSJacob Faibussowitsch     PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, m, N, NULL, &W));
53048e38a8aSPierre Jolivet     /* solve the overlapping 0-block */
5319566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(PC_ApplyOnBlocks, osm->ksp[0], Z, W, 0));
5329566063dSJacob Faibussowitsch     PetscCall(KSPMatSolve(osm->ksp[0], Z, W));
5339566063dSJacob Faibussowitsch     PetscCall(KSPCheckSolve(osm->ksp[0], pc, NULL));
5349566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(PC_ApplyOnBlocks, osm->ksp[0], Z, W, 0));
5359566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&Z));
53648e38a8aSPierre Jolivet 
53748e38a8aSPierre Jolivet     for (i = 0; i < N; ++i) {
5389566063dSJacob Faibussowitsch       PetscCall(VecSet(osm->ly, 0.0));
5399566063dSJacob Faibussowitsch       PetscCall(MatDenseGetColumnVecRead(W, i, &x));
54048e38a8aSPierre Jolivet       if (osm->lprolongation) { /* interpolate the non-overlapping 0-block solution to the local solution (only for restrictive additive) */
5419566063dSJacob Faibussowitsch         PetscCall(VecScatterBegin(osm->lprolongation[0], x, osm->ly, ADD_VALUES, forward));
5429566063dSJacob Faibussowitsch         PetscCall(VecScatterEnd(osm->lprolongation[0], x, osm->ly, ADD_VALUES, forward));
54348e38a8aSPierre Jolivet       } else { /* interpolate the overlapping 0-block solution to the local solution */
5449566063dSJacob Faibussowitsch         PetscCall(VecScatterBegin(osm->lrestriction[0], x, osm->ly, ADD_VALUES, reverse));
5459566063dSJacob Faibussowitsch         PetscCall(VecScatterEnd(osm->lrestriction[0], x, osm->ly, ADD_VALUES, reverse));
54648e38a8aSPierre Jolivet       }
5479566063dSJacob Faibussowitsch       PetscCall(MatDenseRestoreColumnVecRead(W, i, &x));
54848e38a8aSPierre Jolivet 
5499566063dSJacob Faibussowitsch       PetscCall(MatDenseGetColumnVecWrite(Y, i, &x));
55048e38a8aSPierre Jolivet       /* add the local solution to the global solution including the ghost nodes */
5519566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(osm->restriction, osm->ly, x, ADD_VALUES, reverse));
5529566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(osm->restriction, osm->ly, x, ADD_VALUES, reverse));
5539566063dSJacob Faibussowitsch       PetscCall(MatDenseRestoreColumnVecWrite(Y, i, &x));
55448e38a8aSPierre Jolivet     }
5559566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&W));
55698921bdaSJacob Faibussowitsch   } else SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Invalid local composition type: %s", PCCompositeTypes[osm->loctype]);
55748e38a8aSPierre Jolivet   PetscFunctionReturn(0);
55848e38a8aSPierre Jolivet }
55948e38a8aSPierre Jolivet 
5609371c9d4SSatish Balay static PetscErrorCode PCApplyTranspose_ASM(PC pc, Vec x, Vec y) {
5614b9ad928SBarry Smith   PC_ASM     *osm = (PC_ASM *)pc->data;
5621dd8081eSeaulisa   PetscInt    i, n_local_true = osm->n_local_true;
5634b9ad928SBarry Smith   ScatterMode forward = SCATTER_FORWARD, reverse = SCATTER_REVERSE;
5644b9ad928SBarry Smith 
5654b9ad928SBarry Smith   PetscFunctionBegin;
5664b9ad928SBarry Smith   /*
5674b9ad928SBarry Smith      Support for limiting the restriction or interpolation to only local
5684b9ad928SBarry Smith      subdomain values (leaving the other values 0).
5694b9ad928SBarry Smith 
5704b9ad928SBarry Smith      Note: these are reversed from the PCApply_ASM() because we are applying the
5714b9ad928SBarry Smith      transpose of the three terms
5724b9ad928SBarry Smith   */
573d8b3b5e3Seaulisa 
5744b9ad928SBarry Smith   if (!(osm->type & PC_ASM_INTERPOLATE)) {
5754b9ad928SBarry Smith     forward = SCATTER_FORWARD_LOCAL;
5764b9ad928SBarry Smith     /* have to zero the work RHS since scatter may leave some slots empty */
5779566063dSJacob Faibussowitsch     PetscCall(VecSet(osm->lx, 0.0));
5784b9ad928SBarry Smith   }
5792fa5cd67SKarl Rupp   if (!(osm->type & PC_ASM_RESTRICT)) reverse = SCATTER_REVERSE_LOCAL;
5804b9ad928SBarry Smith 
581b0de9f23SBarry Smith   /* zero the global and the local solutions */
5829566063dSJacob Faibussowitsch   PetscCall(VecSet(y, 0.0));
5839566063dSJacob Faibussowitsch   PetscCall(VecSet(osm->ly, 0.0));
584d8b3b5e3Seaulisa 
585b0de9f23SBarry Smith   /* Copy the global RHS to local RHS including the ghost nodes */
5869566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(osm->restriction, x, osm->lx, INSERT_VALUES, forward));
5879566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(osm->restriction, x, osm->lx, INSERT_VALUES, forward));
588d8b3b5e3Seaulisa 
589b0de9f23SBarry Smith   /* Restrict local RHS to the overlapping 0-block RHS */
5909566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward));
5919566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward));
592d8b3b5e3Seaulisa 
5934b9ad928SBarry Smith   /* do the local solves */
594d8b3b5e3Seaulisa   for (i = 0; i < n_local_true; ++i) {
595b0de9f23SBarry Smith     /* solve the overlapping i-block */
5969566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(PC_ApplyOnBlocks, osm->ksp[i], osm->x[i], osm->y[i], 0));
5979566063dSJacob Faibussowitsch     PetscCall(KSPSolveTranspose(osm->ksp[i], osm->x[i], osm->y[i]));
5989566063dSJacob Faibussowitsch     PetscCall(KSPCheckSolve(osm->ksp[i], pc, osm->y[i]));
5999566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(PC_ApplyOnBlocks, osm->ksp[i], osm->x[i], osm->y[i], 0));
600d8b3b5e3Seaulisa 
601910cf402Sprj-     if (osm->lprolongation) { /* interpolate the non-overlapping i-block solution to the local solution */
6029566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward));
6039566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward));
604910cf402Sprj-     } else { /* interpolate the overlapping i-block solution to the local solution */
6059566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse));
6069566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse));
6074b9ad928SBarry Smith     }
608d8b3b5e3Seaulisa 
609d8b3b5e3Seaulisa     if (i < n_local_true - 1) {
610b0de9f23SBarry Smith       /* Restrict local RHS to the overlapping (i+1)-block RHS */
6119566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(osm->lrestriction[i + 1], osm->lx, osm->x[i + 1], INSERT_VALUES, forward));
6129566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(osm->lrestriction[i + 1], osm->lx, osm->x[i + 1], INSERT_VALUES, forward));
6134b9ad928SBarry Smith     }
6144b9ad928SBarry Smith   }
615b0de9f23SBarry Smith   /* Add the local solution to the global solution including the ghost nodes */
6169566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(osm->restriction, osm->ly, y, ADD_VALUES, reverse));
6179566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(osm->restriction, osm->ly, y, ADD_VALUES, reverse));
6184b9ad928SBarry Smith   PetscFunctionReturn(0);
6194b9ad928SBarry Smith }
6204b9ad928SBarry Smith 
6219371c9d4SSatish Balay static PetscErrorCode PCReset_ASM(PC pc) {
6224b9ad928SBarry Smith   PC_ASM  *osm = (PC_ASM *)pc->data;
62313f74950SBarry Smith   PetscInt i;
6244b9ad928SBarry Smith 
6254b9ad928SBarry Smith   PetscFunctionBegin;
62692bb6962SLisandro Dalcin   if (osm->ksp) {
627*48a46eb9SPierre Jolivet     for (i = 0; i < osm->n_local_true; i++) PetscCall(KSPReset(osm->ksp[i]));
62892bb6962SLisandro Dalcin   }
629e09e08ccSBarry Smith   if (osm->pmat) {
630*48a46eb9SPierre Jolivet     if (osm->n_local_true > 0) PetscCall(MatDestroySubMatrices(osm->n_local_true, &osm->pmat));
63192bb6962SLisandro Dalcin   }
6321dd8081eSeaulisa   if (osm->lrestriction) {
6339566063dSJacob Faibussowitsch     PetscCall(VecScatterDestroy(&osm->restriction));
6341dd8081eSeaulisa     for (i = 0; i < osm->n_local_true; i++) {
6359566063dSJacob Faibussowitsch       PetscCall(VecScatterDestroy(&osm->lrestriction[i]));
6369566063dSJacob Faibussowitsch       if (osm->lprolongation) PetscCall(VecScatterDestroy(&osm->lprolongation[i]));
6379566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&osm->x[i]));
6389566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&osm->y[i]));
6394b9ad928SBarry Smith     }
6409566063dSJacob Faibussowitsch     PetscCall(PetscFree(osm->lrestriction));
6419566063dSJacob Faibussowitsch     if (osm->lprolongation) PetscCall(PetscFree(osm->lprolongation));
6429566063dSJacob Faibussowitsch     PetscCall(PetscFree(osm->x));
6439566063dSJacob Faibussowitsch     PetscCall(PetscFree(osm->y));
64492bb6962SLisandro Dalcin   }
6459566063dSJacob Faibussowitsch   PetscCall(PCASMDestroySubdomains(osm->n_local_true, osm->is, osm->is_local));
6469566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&osm->lis));
6479566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&osm->lx));
6489566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&osm->ly));
649*48a46eb9SPierre Jolivet   if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE) PetscCall(MatDestroyMatrices(osm->n_local_true, &osm->lmats));
6502fa5cd67SKarl Rupp 
6519566063dSJacob Faibussowitsch   PetscCall(PetscFree(osm->sub_mat_type));
65280ec0b7dSPatrick Sanan 
6530a545947SLisandro Dalcin   osm->is       = NULL;
6540a545947SLisandro Dalcin   osm->is_local = NULL;
655e91c6855SBarry Smith   PetscFunctionReturn(0);
656e91c6855SBarry Smith }
657e91c6855SBarry Smith 
6589371c9d4SSatish Balay static PetscErrorCode PCDestroy_ASM(PC pc) {
659e91c6855SBarry Smith   PC_ASM  *osm = (PC_ASM *)pc->data;
660e91c6855SBarry Smith   PetscInt i;
661e91c6855SBarry Smith 
662e91c6855SBarry Smith   PetscFunctionBegin;
6639566063dSJacob Faibussowitsch   PetscCall(PCReset_ASM(pc));
664e91c6855SBarry Smith   if (osm->ksp) {
665*48a46eb9SPierre Jolivet     for (i = 0; i < osm->n_local_true; i++) PetscCall(KSPDestroy(&osm->ksp[i]));
6669566063dSJacob Faibussowitsch     PetscCall(PetscFree(osm->ksp));
667e91c6855SBarry Smith   }
6689566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc->data));
66996322394SPierre Jolivet 
6709566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetLocalSubdomains_C", NULL));
6719566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetTotalSubdomains_C", NULL));
6729566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetOverlap_C", NULL));
6739566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetType_C", NULL));
6749566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMGetType_C", NULL));
6759566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetLocalType_C", NULL));
6769566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMGetLocalType_C", NULL));
6779566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetSortIndices_C", NULL));
6789566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMGetSubKSP_C", NULL));
6799566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMGetSubMatType_C", NULL));
6809566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetSubMatType_C", NULL));
6814b9ad928SBarry Smith   PetscFunctionReturn(0);
6824b9ad928SBarry Smith }
6834b9ad928SBarry Smith 
6849371c9d4SSatish Balay static PetscErrorCode PCSetFromOptions_ASM(PC pc, PetscOptionItems *PetscOptionsObject) {
6854b9ad928SBarry Smith   PC_ASM         *osm = (PC_ASM *)pc->data;
6869dcbbd2bSBarry Smith   PetscInt        blocks, ovl;
68757501b6eSBarry Smith   PetscBool       flg;
68892bb6962SLisandro Dalcin   PCASMType       asmtype;
68912cd4985SMatthew G. Knepley   PCCompositeType loctype;
69080ec0b7dSPatrick Sanan   char            sub_mat_type[256];
6914b9ad928SBarry Smith 
6924b9ad928SBarry Smith   PetscFunctionBegin;
693d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Additive Schwarz options");
6949566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_asm_dm_subdomains", "Use DMCreateDomainDecomposition() to define subdomains", "PCASMSetDMSubdomains", osm->dm_subdomains, &osm->dm_subdomains, &flg));
6959566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_asm_blocks", "Number of subdomains", "PCASMSetTotalSubdomains", osm->n, &blocks, &flg));
69665db9045SDmitry Karpeev   if (flg) {
6979566063dSJacob Faibussowitsch     PetscCall(PCASMSetTotalSubdomains(pc, blocks, NULL, NULL));
698d709ea83SDmitry Karpeev     osm->dm_subdomains = PETSC_FALSE;
69965db9045SDmitry Karpeev   }
7009566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_asm_local_blocks", "Number of local subdomains", "PCASMSetLocalSubdomains", osm->n_local_true, &blocks, &flg));
701342c94f9SMatthew G. Knepley   if (flg) {
7029566063dSJacob Faibussowitsch     PetscCall(PCASMSetLocalSubdomains(pc, blocks, NULL, NULL));
703342c94f9SMatthew G. Knepley     osm->dm_subdomains = PETSC_FALSE;
704342c94f9SMatthew G. Knepley   }
7059566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_asm_overlap", "Number of grid points overlap", "PCASMSetOverlap", osm->overlap, &ovl, &flg));
70665db9045SDmitry Karpeev   if (flg) {
7079566063dSJacob Faibussowitsch     PetscCall(PCASMSetOverlap(pc, ovl));
708d709ea83SDmitry Karpeev     osm->dm_subdomains = PETSC_FALSE;
70965db9045SDmitry Karpeev   }
71090d69ab7SBarry Smith   flg = PETSC_FALSE;
7119566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-pc_asm_type", "Type of restriction/extension", "PCASMSetType", PCASMTypes, (PetscEnum)osm->type, (PetscEnum *)&asmtype, &flg));
7129566063dSJacob Faibussowitsch   if (flg) PetscCall(PCASMSetType(pc, asmtype));
71312cd4985SMatthew G. Knepley   flg = PETSC_FALSE;
7149566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-pc_asm_local_type", "Type of local solver composition", "PCASMSetLocalType", PCCompositeTypes, (PetscEnum)osm->loctype, (PetscEnum *)&loctype, &flg));
7159566063dSJacob Faibussowitsch   if (flg) PetscCall(PCASMSetLocalType(pc, loctype));
7169566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFList("-pc_asm_sub_mat_type", "Subsolve Matrix Type", "PCASMSetSubMatType", MatList, NULL, sub_mat_type, 256, &flg));
7171baa6e33SBarry Smith   if (flg) PetscCall(PCASMSetSubMatType(pc, sub_mat_type));
718d0609cedSBarry Smith   PetscOptionsHeadEnd();
7194b9ad928SBarry Smith   PetscFunctionReturn(0);
7204b9ad928SBarry Smith }
7214b9ad928SBarry Smith 
7224b9ad928SBarry Smith /*------------------------------------------------------------------------------------*/
7234b9ad928SBarry Smith 
7249371c9d4SSatish Balay static PetscErrorCode PCASMSetLocalSubdomains_ASM(PC pc, PetscInt n, IS is[], IS is_local[]) {
7254b9ad928SBarry Smith   PC_ASM  *osm = (PC_ASM *)pc->data;
72692bb6962SLisandro Dalcin   PetscInt i;
7274b9ad928SBarry Smith 
7284b9ad928SBarry Smith   PetscFunctionBegin;
72963a3b9bcSJacob Faibussowitsch   PetscCheck(n >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Each process must have 1 or more blocks, n = %" PetscInt_FMT, n);
7307827d75bSBarry Smith   PetscCheck(!pc->setupcalled || (n == osm->n_local_true && !is), PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "PCASMSetLocalSubdomains() should be called before calling PCSetUp().");
731e7e72b3dSBarry Smith 
7324b9ad928SBarry Smith   if (!pc->setupcalled) {
73392bb6962SLisandro Dalcin     if (is) {
7349566063dSJacob Faibussowitsch       for (i = 0; i < n; i++) PetscCall(PetscObjectReference((PetscObject)is[i]));
73592bb6962SLisandro Dalcin     }
736832fc9a5SMatthew Knepley     if (is_local) {
7379566063dSJacob Faibussowitsch       for (i = 0; i < n; i++) PetscCall(PetscObjectReference((PetscObject)is_local[i]));
738832fc9a5SMatthew Knepley     }
7399566063dSJacob Faibussowitsch     PetscCall(PCASMDestroySubdomains(osm->n_local_true, osm->is, osm->is_local));
7402fa5cd67SKarl Rupp 
7414b9ad928SBarry Smith     osm->n_local_true = n;
7420a545947SLisandro Dalcin     osm->is           = NULL;
7430a545947SLisandro Dalcin     osm->is_local     = NULL;
74492bb6962SLisandro Dalcin     if (is) {
7459566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(n, &osm->is));
7462fa5cd67SKarl Rupp       for (i = 0; i < n; i++) osm->is[i] = is[i];
7473d03bb48SJed Brown       /* Flag indicating that the user has set overlapping subdomains so PCASM should not increase their size. */
7483d03bb48SJed Brown       osm->overlap = -1;
74992bb6962SLisandro Dalcin     }
7502b691e39SMatthew Knepley     if (is_local) {
7519566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(n, &osm->is_local));
7522fa5cd67SKarl Rupp       for (i = 0; i < n; i++) osm->is_local[i] = is_local[i];
753a35b7b57SDmitry Karpeev       if (!is) {
7549566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(osm->n_local_true, &osm->is));
755a35b7b57SDmitry Karpeev         for (i = 0; i < osm->n_local_true; i++) {
756a35b7b57SDmitry Karpeev           if (osm->overlap > 0) { /* With positive overlap, osm->is[i] will be modified */
7579566063dSJacob Faibussowitsch             PetscCall(ISDuplicate(osm->is_local[i], &osm->is[i]));
7589566063dSJacob Faibussowitsch             PetscCall(ISCopy(osm->is_local[i], osm->is[i]));
759a35b7b57SDmitry Karpeev           } else {
7609566063dSJacob Faibussowitsch             PetscCall(PetscObjectReference((PetscObject)osm->is_local[i]));
761a35b7b57SDmitry Karpeev             osm->is[i] = osm->is_local[i];
762a35b7b57SDmitry Karpeev           }
763a35b7b57SDmitry Karpeev         }
764a35b7b57SDmitry Karpeev       }
7652b691e39SMatthew Knepley     }
7664b9ad928SBarry Smith   }
7674b9ad928SBarry Smith   PetscFunctionReturn(0);
7684b9ad928SBarry Smith }
7694b9ad928SBarry Smith 
7709371c9d4SSatish Balay static PetscErrorCode PCASMSetTotalSubdomains_ASM(PC pc, PetscInt N, IS *is, IS *is_local) {
7714b9ad928SBarry Smith   PC_ASM     *osm = (PC_ASM *)pc->data;
77213f74950SBarry Smith   PetscMPIInt rank, size;
77378904715SLisandro Dalcin   PetscInt    n;
7744b9ad928SBarry Smith 
7754b9ad928SBarry Smith   PetscFunctionBegin;
77663a3b9bcSJacob Faibussowitsch   PetscCheck(N >= 1, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Number of total blocks must be > 0, N = %" PetscInt_FMT, N);
7777827d75bSBarry Smith   PetscCheck(!is && !is_local, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Use PCASMSetLocalSubdomains() to set specific index sets\n\they cannot be set globally yet.");
7784b9ad928SBarry Smith 
7794b9ad928SBarry Smith   /*
780880770e9SJed Brown      Split the subdomains equally among all processors
7814b9ad928SBarry Smith   */
7829566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank));
7839566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size));
7844b9ad928SBarry Smith   n = N / size + ((N % size) > rank);
78563a3b9bcSJacob Faibussowitsch   PetscCheck(n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Process %d must have at least one block: total processors %d total blocks %" PetscInt_FMT, (int)rank, (int)size, N);
7867827d75bSBarry Smith   PetscCheck(!pc->setupcalled || n == osm->n_local_true, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "PCASMSetTotalSubdomains() should be called before PCSetUp().");
7874b9ad928SBarry Smith   if (!pc->setupcalled) {
7889566063dSJacob Faibussowitsch     PetscCall(PCASMDestroySubdomains(osm->n_local_true, osm->is, osm->is_local));
7892fa5cd67SKarl Rupp 
7904b9ad928SBarry Smith     osm->n_local_true = n;
7910a545947SLisandro Dalcin     osm->is           = NULL;
7920a545947SLisandro Dalcin     osm->is_local     = NULL;
7934b9ad928SBarry Smith   }
7944b9ad928SBarry Smith   PetscFunctionReturn(0);
7954b9ad928SBarry Smith }
7964b9ad928SBarry Smith 
7979371c9d4SSatish Balay static PetscErrorCode PCASMSetOverlap_ASM(PC pc, PetscInt ovl) {
79892bb6962SLisandro Dalcin   PC_ASM *osm = (PC_ASM *)pc->data;
7994b9ad928SBarry Smith 
8004b9ad928SBarry Smith   PetscFunctionBegin;
8017827d75bSBarry Smith   PetscCheck(ovl >= 0, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Negative overlap value requested");
8027827d75bSBarry Smith   PetscCheck(!pc->setupcalled || ovl == osm->overlap, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "PCASMSetOverlap() should be called before PCSetUp().");
8032fa5cd67SKarl Rupp   if (!pc->setupcalled) osm->overlap = ovl;
8044b9ad928SBarry Smith   PetscFunctionReturn(0);
8054b9ad928SBarry Smith }
8064b9ad928SBarry Smith 
8079371c9d4SSatish Balay static PetscErrorCode PCASMSetType_ASM(PC pc, PCASMType type) {
80892bb6962SLisandro Dalcin   PC_ASM *osm = (PC_ASM *)pc->data;
8094b9ad928SBarry Smith 
8104b9ad928SBarry Smith   PetscFunctionBegin;
8114b9ad928SBarry Smith   osm->type     = type;
812bf108f30SBarry Smith   osm->type_set = PETSC_TRUE;
8134b9ad928SBarry Smith   PetscFunctionReturn(0);
8144b9ad928SBarry Smith }
8154b9ad928SBarry Smith 
8169371c9d4SSatish Balay static PetscErrorCode PCASMGetType_ASM(PC pc, PCASMType *type) {
817c60c7ad4SBarry Smith   PC_ASM *osm = (PC_ASM *)pc->data;
818c60c7ad4SBarry Smith 
819c60c7ad4SBarry Smith   PetscFunctionBegin;
820c60c7ad4SBarry Smith   *type = osm->type;
821c60c7ad4SBarry Smith   PetscFunctionReturn(0);
822c60c7ad4SBarry Smith }
823c60c7ad4SBarry Smith 
8249371c9d4SSatish Balay static PetscErrorCode PCASMSetLocalType_ASM(PC pc, PCCompositeType type) {
82512cd4985SMatthew G. Knepley   PC_ASM *osm = (PC_ASM *)pc->data;
82612cd4985SMatthew G. Knepley 
82712cd4985SMatthew G. Knepley   PetscFunctionBegin;
8287827d75bSBarry Smith   PetscCheck(type == PC_COMPOSITE_ADDITIVE || type == PC_COMPOSITE_MULTIPLICATIVE, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Only supports additive or multiplicative as the local type");
82912cd4985SMatthew G. Knepley   osm->loctype = type;
83012cd4985SMatthew G. Knepley   PetscFunctionReturn(0);
83112cd4985SMatthew G. Knepley }
83212cd4985SMatthew G. Knepley 
8339371c9d4SSatish Balay static PetscErrorCode PCASMGetLocalType_ASM(PC pc, PCCompositeType *type) {
83412cd4985SMatthew G. Knepley   PC_ASM *osm = (PC_ASM *)pc->data;
83512cd4985SMatthew G. Knepley 
83612cd4985SMatthew G. Knepley   PetscFunctionBegin;
83712cd4985SMatthew G. Knepley   *type = osm->loctype;
83812cd4985SMatthew G. Knepley   PetscFunctionReturn(0);
83912cd4985SMatthew G. Knepley }
84012cd4985SMatthew G. Knepley 
8419371c9d4SSatish Balay static PetscErrorCode PCASMSetSortIndices_ASM(PC pc, PetscBool doSort) {
8426ed231c7SMatthew Knepley   PC_ASM *osm = (PC_ASM *)pc->data;
8436ed231c7SMatthew Knepley 
8446ed231c7SMatthew Knepley   PetscFunctionBegin;
8456ed231c7SMatthew Knepley   osm->sort_indices = doSort;
8466ed231c7SMatthew Knepley   PetscFunctionReturn(0);
8476ed231c7SMatthew Knepley }
8486ed231c7SMatthew Knepley 
8499371c9d4SSatish Balay static PetscErrorCode PCASMGetSubKSP_ASM(PC pc, PetscInt *n_local, PetscInt *first_local, KSP **ksp) {
85092bb6962SLisandro Dalcin   PC_ASM *osm = (PC_ASM *)pc->data;
8514b9ad928SBarry Smith 
8524b9ad928SBarry Smith   PetscFunctionBegin;
8537827d75bSBarry Smith   PetscCheck(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");
8544b9ad928SBarry Smith 
8552fa5cd67SKarl Rupp   if (n_local) *n_local = osm->n_local_true;
85692bb6962SLisandro Dalcin   if (first_local) {
8579566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Scan(&osm->n_local_true, first_local, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)pc)));
85892bb6962SLisandro Dalcin     *first_local -= osm->n_local_true;
85992bb6962SLisandro Dalcin   }
860ed155784SPierre Jolivet   if (ksp) *ksp = osm->ksp;
8614b9ad928SBarry Smith   PetscFunctionReturn(0);
8624b9ad928SBarry Smith }
8634b9ad928SBarry Smith 
8649371c9d4SSatish Balay static PetscErrorCode PCASMGetSubMatType_ASM(PC pc, MatType *sub_mat_type) {
86580ec0b7dSPatrick Sanan   PC_ASM *osm = (PC_ASM *)pc->data;
86680ec0b7dSPatrick Sanan 
86780ec0b7dSPatrick Sanan   PetscFunctionBegin;
86880ec0b7dSPatrick Sanan   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
86980ec0b7dSPatrick Sanan   PetscValidPointer(sub_mat_type, 2);
87080ec0b7dSPatrick Sanan   *sub_mat_type = osm->sub_mat_type;
87180ec0b7dSPatrick Sanan   PetscFunctionReturn(0);
87280ec0b7dSPatrick Sanan }
87380ec0b7dSPatrick Sanan 
8749371c9d4SSatish Balay static PetscErrorCode PCASMSetSubMatType_ASM(PC pc, MatType sub_mat_type) {
87580ec0b7dSPatrick Sanan   PC_ASM *osm = (PC_ASM *)pc->data;
87680ec0b7dSPatrick Sanan 
87780ec0b7dSPatrick Sanan   PetscFunctionBegin;
87880ec0b7dSPatrick Sanan   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
8799566063dSJacob Faibussowitsch   PetscCall(PetscFree(osm->sub_mat_type));
8809566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(sub_mat_type, (char **)&osm->sub_mat_type));
88180ec0b7dSPatrick Sanan   PetscFunctionReturn(0);
88280ec0b7dSPatrick Sanan }
88380ec0b7dSPatrick Sanan 
8844b9ad928SBarry Smith /*@C
8851093a601SBarry Smith     PCASMSetLocalSubdomains - Sets the local subdomains (for this processor only) for the additive Schwarz preconditioner.
8864b9ad928SBarry Smith 
887d083f849SBarry Smith     Collective on pc
8884b9ad928SBarry Smith 
8894b9ad928SBarry Smith     Input Parameters:
8904b9ad928SBarry Smith +   pc - the preconditioner context
8914b9ad928SBarry Smith .   n - the number of subdomains for this processor (default value = 1)
8928c03b21aSDmitry Karpeev .   is - the index set that defines the subdomains for this processor
8930298fd71SBarry Smith          (or NULL for PETSc to determine subdomains)
894f1ee410cSBarry 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
895f1ee410cSBarry Smith          (or NULL to not provide these)
8964b9ad928SBarry Smith 
897342c94f9SMatthew G. Knepley     Options Database Key:
898342c94f9SMatthew G. Knepley     To set the total number of subdomain blocks rather than specify the
899342c94f9SMatthew G. Knepley     index sets, use the option
900342c94f9SMatthew G. Knepley .    -pc_asm_local_blocks <blks> - Sets local blocks
901342c94f9SMatthew G. Knepley 
9024b9ad928SBarry Smith     Notes:
9031093a601SBarry Smith     The IS numbering is in the parallel, global numbering of the vector for both is and is_local
9044b9ad928SBarry Smith 
9054b9ad928SBarry Smith     By default the ASM preconditioner uses 1 block per processor.
9064b9ad928SBarry Smith 
9074b9ad928SBarry Smith     Use PCASMSetTotalSubdomains() to set the subdomains for all processors.
9084b9ad928SBarry Smith 
909f1ee410cSBarry Smith     If is_local is provided and PCASMType is PC_ASM_RESTRICT then the solution only over the is_local region is interpolated
910f1ee410cSBarry Smith     back to form the global solution (this is the standard restricted additive Schwarz method)
911f1ee410cSBarry Smith 
912f1ee410cSBarry 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
913f1ee410cSBarry Smith     no code to handle that case.
914f1ee410cSBarry Smith 
9154b9ad928SBarry Smith     Level: advanced
9164b9ad928SBarry Smith 
917db781477SPatrick Sanan .seealso: `PCASMSetTotalSubdomains()`, `PCASMSetOverlap()`, `PCASMGetSubKSP()`,
918db781477SPatrick Sanan           `PCASMCreateSubdomains2D()`, `PCASMGetLocalSubdomains()`, `PCASMType`, `PCASMSetType()`
9194b9ad928SBarry Smith @*/
9209371c9d4SSatish Balay PetscErrorCode PCASMSetLocalSubdomains(PC pc, PetscInt n, IS is[], IS is_local[]) {
9214b9ad928SBarry Smith   PetscFunctionBegin;
9220700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
923cac4c232SBarry Smith   PetscTryMethod(pc, "PCASMSetLocalSubdomains_C", (PC, PetscInt, IS[], IS[]), (pc, n, is, is_local));
9244b9ad928SBarry Smith   PetscFunctionReturn(0);
9254b9ad928SBarry Smith }
9264b9ad928SBarry Smith 
9274b9ad928SBarry Smith /*@C
928feb221f8SDmitry Karpeev     PCASMSetTotalSubdomains - Sets the subdomains for all processors for the
9294b9ad928SBarry Smith     additive Schwarz preconditioner.  Either all or no processors in the
9304b9ad928SBarry Smith     PC communicator must call this routine, with the same index sets.
9314b9ad928SBarry Smith 
932d083f849SBarry Smith     Collective on pc
9334b9ad928SBarry Smith 
9344b9ad928SBarry Smith     Input Parameters:
9354b9ad928SBarry Smith +   pc - the preconditioner context
936feb221f8SDmitry Karpeev .   N  - the number of subdomains for all processors
937feb221f8SDmitry Karpeev .   is - the index sets that define the subdomains for all processors
938dfaa0c5fSBarry Smith          (or NULL to ask PETSc to determine the subdomains)
9392b691e39SMatthew Knepley -   is_local - the index sets that define the local part of the subdomains for this processor
940f1ee410cSBarry Smith          (or NULL to not provide this information)
9414b9ad928SBarry Smith 
9424b9ad928SBarry Smith     Options Database Key:
9434b9ad928SBarry Smith     To set the total number of subdomain blocks rather than specify the
9444b9ad928SBarry Smith     index sets, use the option
9454b9ad928SBarry Smith .    -pc_asm_blocks <blks> - Sets total blocks
9464b9ad928SBarry Smith 
9474b9ad928SBarry Smith     Notes:
948f1ee410cSBarry Smith     Currently you cannot use this to set the actual subdomains with the argument is or is_local.
9494b9ad928SBarry Smith 
9504b9ad928SBarry Smith     By default the ASM preconditioner uses 1 block per processor.
9514b9ad928SBarry Smith 
9524b9ad928SBarry Smith     These index sets cannot be destroyed until after completion of the
9534b9ad928SBarry Smith     linear solves for which the ASM preconditioner is being used.
9544b9ad928SBarry Smith 
9554b9ad928SBarry Smith     Use PCASMSetLocalSubdomains() to set local subdomains.
9564b9ad928SBarry Smith 
9571093a601SBarry Smith     The IS numbering is in the parallel, global numbering of the vector for both is and is_local
9581093a601SBarry Smith 
9594b9ad928SBarry Smith     Level: advanced
9604b9ad928SBarry Smith 
961db781477SPatrick Sanan .seealso: `PCASMSetLocalSubdomains()`, `PCASMSetOverlap()`, `PCASMGetSubKSP()`,
962db781477SPatrick Sanan           `PCASMCreateSubdomains2D()`
9634b9ad928SBarry Smith @*/
9649371c9d4SSatish Balay PetscErrorCode PCASMSetTotalSubdomains(PC pc, PetscInt N, IS is[], IS is_local[]) {
9654b9ad928SBarry Smith   PetscFunctionBegin;
9660700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
967cac4c232SBarry Smith   PetscTryMethod(pc, "PCASMSetTotalSubdomains_C", (PC, PetscInt, IS[], IS[]), (pc, N, is, is_local));
9684b9ad928SBarry Smith   PetscFunctionReturn(0);
9694b9ad928SBarry Smith }
9704b9ad928SBarry Smith 
9714b9ad928SBarry Smith /*@
9724b9ad928SBarry Smith     PCASMSetOverlap - Sets the overlap between a pair of subdomains for the
9734b9ad928SBarry Smith     additive Schwarz preconditioner.  Either all or no processors in the
974f1ee410cSBarry Smith     PC communicator must call this routine.
9754b9ad928SBarry Smith 
976d083f849SBarry Smith     Logically Collective on pc
9774b9ad928SBarry Smith 
9784b9ad928SBarry Smith     Input Parameters:
9794b9ad928SBarry Smith +   pc  - the preconditioner context
9804b9ad928SBarry Smith -   ovl - the amount of overlap between subdomains (ovl >= 0, default value = 1)
9814b9ad928SBarry Smith 
9824b9ad928SBarry Smith     Options Database Key:
9834b9ad928SBarry Smith .   -pc_asm_overlap <ovl> - Sets overlap
9844b9ad928SBarry Smith 
9854b9ad928SBarry Smith     Notes:
9864b9ad928SBarry Smith     By default the ASM preconditioner uses 1 block per processor.  To use
9874b9ad928SBarry Smith     multiple blocks per perocessor, see PCASMSetTotalSubdomains() and
9884b9ad928SBarry Smith     PCASMSetLocalSubdomains() (and the option -pc_asm_blocks <blks>).
9894b9ad928SBarry Smith 
9904b9ad928SBarry Smith     The overlap defaults to 1, so if one desires that no additional
9914b9ad928SBarry Smith     overlap be computed beyond what may have been set with a call to
9924b9ad928SBarry Smith     PCASMSetTotalSubdomains() or PCASMSetLocalSubdomains(), then ovl
9934b9ad928SBarry Smith     must be set to be 0.  In particular, if one does not explicitly set
9944b9ad928SBarry Smith     the subdomains an application code, then all overlap would be computed
9954b9ad928SBarry Smith     internally by PETSc, and using an overlap of 0 would result in an ASM
9964b9ad928SBarry Smith     variant that is equivalent to the block Jacobi preconditioner.
9974b9ad928SBarry Smith 
998f1ee410cSBarry Smith     The default algorithm used by PETSc to increase overlap is fast, but not scalable,
999f1ee410cSBarry Smith     use the option -mat_increase_overlap_scalable when the problem and number of processes is large.
1000f1ee410cSBarry Smith 
10014b9ad928SBarry Smith     Note that one can define initial index sets with any overlap via
1002f1ee410cSBarry Smith     PCASMSetLocalSubdomains(); the routine
10034b9ad928SBarry Smith     PCASMSetOverlap() merely allows PETSc to extend that overlap further
10044b9ad928SBarry Smith     if desired.
10054b9ad928SBarry Smith 
10064b9ad928SBarry Smith     Level: intermediate
10074b9ad928SBarry Smith 
1008db781477SPatrick Sanan .seealso: `PCASMSetTotalSubdomains()`, `PCASMSetLocalSubdomains()`, `PCASMGetSubKSP()`,
1009db781477SPatrick Sanan           `PCASMCreateSubdomains2D()`, `PCASMGetLocalSubdomains()`, `MatIncreaseOverlap()`
10104b9ad928SBarry Smith @*/
10119371c9d4SSatish Balay PetscErrorCode PCASMSetOverlap(PC pc, PetscInt ovl) {
10124b9ad928SBarry Smith   PetscFunctionBegin;
10130700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1014c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc, ovl, 2);
1015cac4c232SBarry Smith   PetscTryMethod(pc, "PCASMSetOverlap_C", (PC, PetscInt), (pc, ovl));
10164b9ad928SBarry Smith   PetscFunctionReturn(0);
10174b9ad928SBarry Smith }
10184b9ad928SBarry Smith 
10194b9ad928SBarry Smith /*@
10204b9ad928SBarry Smith     PCASMSetType - Sets the type of restriction and interpolation used
10214b9ad928SBarry Smith     for local problems in the additive Schwarz method.
10224b9ad928SBarry Smith 
1023d083f849SBarry Smith     Logically Collective on pc
10244b9ad928SBarry Smith 
10254b9ad928SBarry Smith     Input Parameters:
10264b9ad928SBarry Smith +   pc  - the preconditioner context
10274b9ad928SBarry Smith -   type - variant of ASM, one of
10284b9ad928SBarry Smith .vb
10294b9ad928SBarry Smith       PC_ASM_BASIC       - full interpolation and restriction
103082b5ce2aSStefano Zampini       PC_ASM_RESTRICT    - full restriction, local processor interpolation (default)
10314b9ad928SBarry Smith       PC_ASM_INTERPOLATE - full interpolation, local processor restriction
10324b9ad928SBarry Smith       PC_ASM_NONE        - local processor restriction and interpolation
10334b9ad928SBarry Smith .ve
10344b9ad928SBarry Smith 
10354b9ad928SBarry Smith     Options Database Key:
10364b9ad928SBarry Smith .   -pc_asm_type [basic,restrict,interpolate,none] - Sets ASM type
10374b9ad928SBarry Smith 
103895452b02SPatrick Sanan     Notes:
103995452b02SPatrick Sanan     if the is_local arguments are passed to PCASMSetLocalSubdomains() then they are used when PC_ASM_RESTRICT has been selected
1040f1ee410cSBarry Smith     to limit the local processor interpolation
1041f1ee410cSBarry Smith 
10424b9ad928SBarry Smith     Level: intermediate
10434b9ad928SBarry Smith 
1044db781477SPatrick Sanan .seealso: `PCASMSetTotalSubdomains()`, `PCASMSetTotalSubdomains()`, `PCASMGetSubKSP()`,
1045db781477SPatrick Sanan           `PCASMCreateSubdomains2D()`, `PCASMType`, `PCASMSetLocalType()`, `PCASMGetLocalType()`
10464b9ad928SBarry Smith @*/
10479371c9d4SSatish Balay PetscErrorCode PCASMSetType(PC pc, PCASMType type) {
10484b9ad928SBarry Smith   PetscFunctionBegin;
10490700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1050c5eb9154SBarry Smith   PetscValidLogicalCollectiveEnum(pc, type, 2);
1051cac4c232SBarry Smith   PetscTryMethod(pc, "PCASMSetType_C", (PC, PCASMType), (pc, type));
10524b9ad928SBarry Smith   PetscFunctionReturn(0);
10534b9ad928SBarry Smith }
10544b9ad928SBarry Smith 
1055c60c7ad4SBarry Smith /*@
1056c60c7ad4SBarry Smith     PCASMGetType - Gets the type of restriction and interpolation used
1057c60c7ad4SBarry Smith     for local problems in the additive Schwarz method.
1058c60c7ad4SBarry Smith 
1059d083f849SBarry Smith     Logically Collective on pc
1060c60c7ad4SBarry Smith 
1061c60c7ad4SBarry Smith     Input Parameter:
1062c60c7ad4SBarry Smith .   pc  - the preconditioner context
1063c60c7ad4SBarry Smith 
1064c60c7ad4SBarry Smith     Output Parameter:
1065c60c7ad4SBarry Smith .   type - variant of ASM, one of
1066c60c7ad4SBarry Smith 
1067c60c7ad4SBarry Smith .vb
1068c60c7ad4SBarry Smith       PC_ASM_BASIC       - full interpolation and restriction
1069c60c7ad4SBarry Smith       PC_ASM_RESTRICT    - full restriction, local processor interpolation
1070c60c7ad4SBarry Smith       PC_ASM_INTERPOLATE - full interpolation, local processor restriction
1071c60c7ad4SBarry Smith       PC_ASM_NONE        - local processor restriction and interpolation
1072c60c7ad4SBarry Smith .ve
1073c60c7ad4SBarry Smith 
1074c60c7ad4SBarry Smith     Options Database Key:
1075c60c7ad4SBarry Smith .   -pc_asm_type [basic,restrict,interpolate,none] - Sets ASM type
1076c60c7ad4SBarry Smith 
1077c60c7ad4SBarry Smith     Level: intermediate
1078c60c7ad4SBarry Smith 
1079db781477SPatrick Sanan .seealso: `PCASMSetTotalSubdomains()`, `PCASMSetTotalSubdomains()`, `PCASMGetSubKSP()`,
1080db781477SPatrick Sanan           `PCASMCreateSubdomains2D()`, `PCASMType`, `PCASMSetType()`, `PCASMSetLocalType()`, `PCASMGetLocalType()`
1081c60c7ad4SBarry Smith @*/
10829371c9d4SSatish Balay PetscErrorCode PCASMGetType(PC pc, PCASMType *type) {
1083c60c7ad4SBarry Smith   PetscFunctionBegin;
1084c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1085cac4c232SBarry Smith   PetscUseMethod(pc, "PCASMGetType_C", (PC, PCASMType *), (pc, type));
1086c60c7ad4SBarry Smith   PetscFunctionReturn(0);
1087c60c7ad4SBarry Smith }
1088c60c7ad4SBarry Smith 
108912cd4985SMatthew G. Knepley /*@
109012cd4985SMatthew G. Knepley   PCASMSetLocalType - Sets the type of composition used for local problems in the additive Schwarz method.
109112cd4985SMatthew G. Knepley 
1092d083f849SBarry Smith   Logically Collective on pc
109312cd4985SMatthew G. Knepley 
109412cd4985SMatthew G. Knepley   Input Parameters:
109512cd4985SMatthew G. Knepley + pc  - the preconditioner context
109612cd4985SMatthew G. Knepley - type - type of composition, one of
109712cd4985SMatthew G. Knepley .vb
109812cd4985SMatthew G. Knepley   PC_COMPOSITE_ADDITIVE       - local additive combination
109912cd4985SMatthew G. Knepley   PC_COMPOSITE_MULTIPLICATIVE - local multiplicative combination
110012cd4985SMatthew G. Knepley .ve
110112cd4985SMatthew G. Knepley 
110212cd4985SMatthew G. Knepley   Options Database Key:
110312cd4985SMatthew G. Knepley . -pc_asm_local_type [additive,multiplicative] - Sets local solver composition type
110412cd4985SMatthew G. Knepley 
110512cd4985SMatthew G. Knepley   Level: intermediate
110612cd4985SMatthew G. Knepley 
1107db781477SPatrick Sanan .seealso: `PCASMSetType()`, `PCASMGetType()`, `PCASMGetLocalType()`, `PCASM`, `PCASMType`, `PCASMSetType()`, `PCASMGetType()`, `PCCompositeType`
110812cd4985SMatthew G. Knepley @*/
11099371c9d4SSatish Balay PetscErrorCode PCASMSetLocalType(PC pc, PCCompositeType type) {
111012cd4985SMatthew G. Knepley   PetscFunctionBegin;
111112cd4985SMatthew G. Knepley   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
111212cd4985SMatthew G. Knepley   PetscValidLogicalCollectiveEnum(pc, type, 2);
1113cac4c232SBarry Smith   PetscTryMethod(pc, "PCASMSetLocalType_C", (PC, PCCompositeType), (pc, type));
111412cd4985SMatthew G. Knepley   PetscFunctionReturn(0);
111512cd4985SMatthew G. Knepley }
111612cd4985SMatthew G. Knepley 
111712cd4985SMatthew G. Knepley /*@
111812cd4985SMatthew G. Knepley   PCASMGetLocalType - Gets the type of composition used for local problems in the additive Schwarz method.
111912cd4985SMatthew G. Knepley 
1120d083f849SBarry Smith   Logically Collective on pc
112112cd4985SMatthew G. Knepley 
112212cd4985SMatthew G. Knepley   Input Parameter:
112312cd4985SMatthew G. Knepley . pc  - the preconditioner context
112412cd4985SMatthew G. Knepley 
112512cd4985SMatthew G. Knepley   Output Parameter:
112612cd4985SMatthew G. Knepley . type - type of composition, one of
112712cd4985SMatthew G. Knepley .vb
112812cd4985SMatthew G. Knepley   PC_COMPOSITE_ADDITIVE       - local additive combination
112912cd4985SMatthew G. Knepley   PC_COMPOSITE_MULTIPLICATIVE - local multiplicative combination
113012cd4985SMatthew G. Knepley .ve
113112cd4985SMatthew G. Knepley 
113212cd4985SMatthew G. Knepley   Options Database Key:
113312cd4985SMatthew G. Knepley . -pc_asm_local_type [additive,multiplicative] - Sets local solver composition type
113412cd4985SMatthew G. Knepley 
113512cd4985SMatthew G. Knepley   Level: intermediate
113612cd4985SMatthew G. Knepley 
1137db781477SPatrick Sanan .seealso: `PCASMSetType()`, `PCASMGetType()`, `PCASMSetLocalType()`, `PCASMCreate()`, `PCASMType`, `PCASMSetType()`, `PCASMGetType()`, `PCCompositeType`
113812cd4985SMatthew G. Knepley @*/
11399371c9d4SSatish Balay PetscErrorCode PCASMGetLocalType(PC pc, PCCompositeType *type) {
114012cd4985SMatthew G. Knepley   PetscFunctionBegin;
114112cd4985SMatthew G. Knepley   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
114212cd4985SMatthew G. Knepley   PetscValidPointer(type, 2);
1143cac4c232SBarry Smith   PetscUseMethod(pc, "PCASMGetLocalType_C", (PC, PCCompositeType *), (pc, type));
114412cd4985SMatthew G. Knepley   PetscFunctionReturn(0);
114512cd4985SMatthew G. Knepley }
114612cd4985SMatthew G. Knepley 
11476ed231c7SMatthew Knepley /*@
11486ed231c7SMatthew Knepley     PCASMSetSortIndices - Determines whether subdomain indices are sorted.
11496ed231c7SMatthew Knepley 
1150d083f849SBarry Smith     Logically Collective on pc
11516ed231c7SMatthew Knepley 
11526ed231c7SMatthew Knepley     Input Parameters:
11536ed231c7SMatthew Knepley +   pc  - the preconditioner context
11546ed231c7SMatthew Knepley -   doSort - sort the subdomain indices
11556ed231c7SMatthew Knepley 
11566ed231c7SMatthew Knepley     Level: intermediate
11576ed231c7SMatthew Knepley 
1158db781477SPatrick Sanan .seealso: `PCASMSetLocalSubdomains()`, `PCASMSetTotalSubdomains()`, `PCASMGetSubKSP()`,
1159db781477SPatrick Sanan           `PCASMCreateSubdomains2D()`
11606ed231c7SMatthew Knepley @*/
11619371c9d4SSatish Balay PetscErrorCode PCASMSetSortIndices(PC pc, PetscBool doSort) {
11626ed231c7SMatthew Knepley   PetscFunctionBegin;
11630700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1164acfcf0e5SJed Brown   PetscValidLogicalCollectiveBool(pc, doSort, 2);
1165cac4c232SBarry Smith   PetscTryMethod(pc, "PCASMSetSortIndices_C", (PC, PetscBool), (pc, doSort));
11666ed231c7SMatthew Knepley   PetscFunctionReturn(0);
11676ed231c7SMatthew Knepley }
11686ed231c7SMatthew Knepley 
11694b9ad928SBarry Smith /*@C
11704b9ad928SBarry Smith    PCASMGetSubKSP - Gets the local KSP contexts for all blocks on
11714b9ad928SBarry Smith    this processor.
11724b9ad928SBarry Smith 
1173d083f849SBarry Smith    Collective on pc iff first_local is requested
11744b9ad928SBarry Smith 
11754b9ad928SBarry Smith    Input Parameter:
11764b9ad928SBarry Smith .  pc - the preconditioner context
11774b9ad928SBarry Smith 
11784b9ad928SBarry Smith    Output Parameters:
11790298fd71SBarry Smith +  n_local - the number of blocks on this processor or NULL
11800298fd71SBarry Smith .  first_local - the global number of the first block on this processor or NULL,
11810298fd71SBarry Smith                  all processors must request or all must pass NULL
11824b9ad928SBarry Smith -  ksp - the array of KSP contexts
11834b9ad928SBarry Smith 
11844b9ad928SBarry Smith    Note:
1185d29017ddSJed Brown    After PCASMGetSubKSP() the array of KSPes is not to be freed.
11864b9ad928SBarry Smith 
11874b9ad928SBarry Smith    You must call KSPSetUp() before calling PCASMGetSubKSP().
11884b9ad928SBarry Smith 
1189d29017ddSJed Brown    Fortran note:
11902bf68e3eSBarry 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.
1191d29017ddSJed Brown 
11924b9ad928SBarry Smith    Level: advanced
11934b9ad928SBarry Smith 
1194db781477SPatrick Sanan .seealso: `PCASMSetTotalSubdomains()`, `PCASMSetTotalSubdomains()`, `PCASMSetOverlap()`,
1195db781477SPatrick Sanan           `PCASMCreateSubdomains2D()`,
11964b9ad928SBarry Smith @*/
11979371c9d4SSatish Balay PetscErrorCode PCASMGetSubKSP(PC pc, PetscInt *n_local, PetscInt *first_local, KSP *ksp[]) {
11984b9ad928SBarry Smith   PetscFunctionBegin;
11990700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1200cac4c232SBarry Smith   PetscUseMethod(pc, "PCASMGetSubKSP_C", (PC, PetscInt *, PetscInt *, KSP **), (pc, n_local, first_local, ksp));
12014b9ad928SBarry Smith   PetscFunctionReturn(0);
12024b9ad928SBarry Smith }
12034b9ad928SBarry Smith 
12044b9ad928SBarry Smith /* -------------------------------------------------------------------------------------*/
12054b9ad928SBarry Smith /*MC
12063b09bd56SBarry Smith    PCASM - Use the (restricted) additive Schwarz method, each block is (approximately) solved with
12074b9ad928SBarry Smith            its own KSP object.
12084b9ad928SBarry Smith 
12094b9ad928SBarry Smith    Options Database Keys:
121049517cdeSBarry Smith +  -pc_asm_blocks <blks> - Sets total blocks
12114b9ad928SBarry Smith .  -pc_asm_overlap <ovl> - Sets overlap
1212f1ee410cSBarry Smith .  -pc_asm_type [basic,restrict,interpolate,none] - Sets ASM type, default is restrict
1213f1ee410cSBarry Smith -  -pc_asm_local_type [additive, multiplicative] - Sets ASM type, default is additive
12144b9ad928SBarry Smith 
12153b09bd56SBarry Smith      IMPORTANT: If you run with, for example, 3 blocks on 1 processor or 3 blocks on 3 processors you
12163b09bd56SBarry Smith       will get a different convergence rate due to the default option of -pc_asm_type restrict. Use
12173b09bd56SBarry Smith       -pc_asm_type basic to use the standard ASM.
12183b09bd56SBarry Smith 
121995452b02SPatrick Sanan    Notes:
122095452b02SPatrick Sanan     Each processor can have one or more blocks, but a block cannot be shared by more
1221f1ee410cSBarry Smith      than one processor. Use PCGASM for subdomains shared by multiple processes. Defaults to one block per processor.
12224b9ad928SBarry Smith 
12233b09bd56SBarry Smith      To set options on the solvers for each block append -sub_ to all the KSP, and PC
1224d7ee0231SBarry Smith         options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly
12254b9ad928SBarry Smith 
1226a8c7a070SBarry Smith      To set the options on the solvers separate for each block call PCASMGetSubKSP()
12274b9ad928SBarry Smith          and set the options directly on the resulting KSP object (you can access its PC
12284b9ad928SBarry Smith          with KSPGetPC())
12294b9ad928SBarry Smith 
12304b9ad928SBarry Smith    Level: beginner
12314b9ad928SBarry Smith 
1232c582cd25SBarry Smith     References:
1233606c0280SSatish Balay +   * - M Dryja, OB Widlund, An additive variant of the Schwarz alternating method for the case of many subregions
123496a0c994SBarry Smith      Courant Institute, New York University Technical report
1235606c0280SSatish Balay -   * - Barry Smith, Petter Bjorstad, and William Gropp, Domain Decompositions: Parallel Multilevel Methods for Elliptic Partial Differential Equations,
123696a0c994SBarry Smith     Cambridge University Press.
1237c582cd25SBarry Smith 
1238db781477SPatrick Sanan .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`,
1239db781477SPatrick Sanan           `PCBJACOBI`, `PCASMGetSubKSP()`, `PCASMSetLocalSubdomains()`, `PCASMType`, `PCASMGetType()`, `PCASMSetLocalType()`, `PCASMGetLocalType()`
1240db781477SPatrick Sanan           `PCASMSetTotalSubdomains()`, `PCSetModifySubMatrices()`, `PCASMSetOverlap()`, `PCASMSetType()`, `PCCompositeType`
1241e09e08ccSBarry Smith 
12424b9ad928SBarry Smith M*/
12434b9ad928SBarry Smith 
12449371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode PCCreate_ASM(PC pc) {
12454b9ad928SBarry Smith   PC_ASM *osm;
12464b9ad928SBarry Smith 
12474b9ad928SBarry Smith   PetscFunctionBegin;
12489566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(pc, &osm));
12492fa5cd67SKarl Rupp 
12504b9ad928SBarry Smith   osm->n             = PETSC_DECIDE;
12514b9ad928SBarry Smith   osm->n_local       = 0;
12522b837212SDmitry Karpeev   osm->n_local_true  = PETSC_DECIDE;
12534b9ad928SBarry Smith   osm->overlap       = 1;
12540a545947SLisandro Dalcin   osm->ksp           = NULL;
12550a545947SLisandro Dalcin   osm->restriction   = NULL;
12560a545947SLisandro Dalcin   osm->lprolongation = NULL;
12570a545947SLisandro Dalcin   osm->lrestriction  = NULL;
12580a545947SLisandro Dalcin   osm->x             = NULL;
12590a545947SLisandro Dalcin   osm->y             = NULL;
12600a545947SLisandro Dalcin   osm->is            = NULL;
12610a545947SLisandro Dalcin   osm->is_local      = NULL;
12620a545947SLisandro Dalcin   osm->mat           = NULL;
12630a545947SLisandro Dalcin   osm->pmat          = NULL;
12644b9ad928SBarry Smith   osm->type          = PC_ASM_RESTRICT;
126512cd4985SMatthew G. Knepley   osm->loctype       = PC_COMPOSITE_ADDITIVE;
12666ed231c7SMatthew Knepley   osm->sort_indices  = PETSC_TRUE;
1267d709ea83SDmitry Karpeev   osm->dm_subdomains = PETSC_FALSE;
126880ec0b7dSPatrick Sanan   osm->sub_mat_type  = NULL;
12694b9ad928SBarry Smith 
127092bb6962SLisandro Dalcin   pc->data                 = (void *)osm;
12714b9ad928SBarry Smith   pc->ops->apply           = PCApply_ASM;
127248e38a8aSPierre Jolivet   pc->ops->matapply        = PCMatApply_ASM;
12734b9ad928SBarry Smith   pc->ops->applytranspose  = PCApplyTranspose_ASM;
12744b9ad928SBarry Smith   pc->ops->setup           = PCSetUp_ASM;
1275e91c6855SBarry Smith   pc->ops->reset           = PCReset_ASM;
12764b9ad928SBarry Smith   pc->ops->destroy         = PCDestroy_ASM;
12774b9ad928SBarry Smith   pc->ops->setfromoptions  = PCSetFromOptions_ASM;
12784b9ad928SBarry Smith   pc->ops->setuponblocks   = PCSetUpOnBlocks_ASM;
12794b9ad928SBarry Smith   pc->ops->view            = PCView_ASM;
12800a545947SLisandro Dalcin   pc->ops->applyrichardson = NULL;
12814b9ad928SBarry Smith 
12829566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetLocalSubdomains_C", PCASMSetLocalSubdomains_ASM));
12839566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetTotalSubdomains_C", PCASMSetTotalSubdomains_ASM));
12849566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetOverlap_C", PCASMSetOverlap_ASM));
12859566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetType_C", PCASMSetType_ASM));
12869566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMGetType_C", PCASMGetType_ASM));
12879566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetLocalType_C", PCASMSetLocalType_ASM));
12889566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMGetLocalType_C", PCASMGetLocalType_ASM));
12899566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetSortIndices_C", PCASMSetSortIndices_ASM));
12909566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMGetSubKSP_C", PCASMGetSubKSP_ASM));
12919566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMGetSubMatType_C", PCASMGetSubMatType_ASM));
12929566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetSubMatType_C", PCASMSetSubMatType_ASM));
12934b9ad928SBarry Smith   PetscFunctionReturn(0);
12944b9ad928SBarry Smith }
12954b9ad928SBarry Smith 
129692bb6962SLisandro Dalcin /*@C
129792bb6962SLisandro Dalcin    PCASMCreateSubdomains - Creates the index sets for the overlapping Schwarz
12983b8d980fSPierre Jolivet    preconditioner for any problem on a general grid.
129992bb6962SLisandro Dalcin 
130092bb6962SLisandro Dalcin    Collective
130192bb6962SLisandro Dalcin 
130292bb6962SLisandro Dalcin    Input Parameters:
130392bb6962SLisandro Dalcin +  A - The global matrix operator
130492bb6962SLisandro Dalcin -  n - the number of local blocks
130592bb6962SLisandro Dalcin 
130692bb6962SLisandro Dalcin    Output Parameters:
130792bb6962SLisandro Dalcin .  outis - the array of index sets defining the subdomains
130892bb6962SLisandro Dalcin 
130992bb6962SLisandro Dalcin    Level: advanced
131092bb6962SLisandro Dalcin 
13117d6bfa3bSBarry Smith    Note: this generates nonoverlapping subdomains; the PCASM will generate the overlap
13127d6bfa3bSBarry Smith     from these if you use PCASMSetLocalSubdomains()
13137d6bfa3bSBarry Smith 
13147d6bfa3bSBarry Smith     In the Fortran version you must provide the array outis[] already allocated of length n.
13157d6bfa3bSBarry Smith 
1316db781477SPatrick Sanan .seealso: `PCASMSetLocalSubdomains()`, `PCASMDestroySubdomains()`
131792bb6962SLisandro Dalcin @*/
13189371c9d4SSatish Balay PetscErrorCode PCASMCreateSubdomains(Mat A, PetscInt n, IS *outis[]) {
131992bb6962SLisandro Dalcin   MatPartitioning mpart;
132092bb6962SLisandro Dalcin   const char     *prefix;
132192bb6962SLisandro Dalcin   PetscInt        i, j, rstart, rend, bs;
1322976e8c5aSLisandro Dalcin   PetscBool       hasop, isbaij = PETSC_FALSE, foundpart = PETSC_FALSE;
13230298fd71SBarry Smith   Mat             Ad = NULL, adj;
132492bb6962SLisandro Dalcin   IS              ispart, isnumb, *is;
132592bb6962SLisandro Dalcin 
132692bb6962SLisandro Dalcin   PetscFunctionBegin;
13270700a824SBarry Smith   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
132892bb6962SLisandro Dalcin   PetscValidPointer(outis, 3);
132963a3b9bcSJacob Faibussowitsch   PetscCheck(n >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "number of local blocks must be > 0, n = %" PetscInt_FMT, n);
133092bb6962SLisandro Dalcin 
133192bb6962SLisandro Dalcin   /* Get prefix, row distribution, and block size */
13329566063dSJacob Faibussowitsch   PetscCall(MatGetOptionsPrefix(A, &prefix));
13339566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(A, &rstart, &rend));
13349566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(A, &bs));
133563a3b9bcSJacob Faibussowitsch   PetscCheck(rstart / bs * bs == rstart && rend / bs * bs == rend, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "bad row distribution [%" PetscInt_FMT ",%" PetscInt_FMT ") for matrix block size %" PetscInt_FMT, rstart, rend, bs);
133665e19b50SBarry Smith 
133792bb6962SLisandro Dalcin   /* Get diagonal block from matrix if possible */
13389566063dSJacob Faibussowitsch   PetscCall(MatHasOperation(A, MATOP_GET_DIAGONAL_BLOCK, &hasop));
1339*48a46eb9SPierre Jolivet   if (hasop) PetscCall(MatGetDiagonalBlock(A, &Ad));
134092bb6962SLisandro Dalcin   if (Ad) {
13419566063dSJacob Faibussowitsch     PetscCall(PetscObjectBaseTypeCompare((PetscObject)Ad, MATSEQBAIJ, &isbaij));
13429566063dSJacob Faibussowitsch     if (!isbaij) PetscCall(PetscObjectBaseTypeCompare((PetscObject)Ad, MATSEQSBAIJ, &isbaij));
134392bb6962SLisandro Dalcin   }
134492bb6962SLisandro Dalcin   if (Ad && n > 1) {
1345ace3abfcSBarry Smith     PetscBool match, done;
134692bb6962SLisandro Dalcin     /* Try to setup a good matrix partitioning if available */
13479566063dSJacob Faibussowitsch     PetscCall(MatPartitioningCreate(PETSC_COMM_SELF, &mpart));
13489566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)mpart, prefix));
13499566063dSJacob Faibussowitsch     PetscCall(MatPartitioningSetFromOptions(mpart));
13509566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)mpart, MATPARTITIONINGCURRENT, &match));
1351*48a46eb9SPierre Jolivet     if (!match) PetscCall(PetscObjectTypeCompare((PetscObject)mpart, MATPARTITIONINGSQUARE, &match));
135292bb6962SLisandro Dalcin     if (!match) { /* assume a "good" partitioner is available */
13531a83f524SJed Brown       PetscInt        na;
13541a83f524SJed Brown       const PetscInt *ia, *ja;
13559566063dSJacob Faibussowitsch       PetscCall(MatGetRowIJ(Ad, 0, PETSC_TRUE, isbaij, &na, &ia, &ja, &done));
135692bb6962SLisandro Dalcin       if (done) {
135792bb6962SLisandro Dalcin         /* Build adjacency matrix by hand. Unfortunately a call to
135892bb6962SLisandro Dalcin            MatConvert(Ad,MATMPIADJ,MAT_INITIAL_MATRIX,&adj) will
135992bb6962SLisandro Dalcin            remove the block-aij structure and we cannot expect
136092bb6962SLisandro Dalcin            MatPartitioning to split vertices as we need */
13610a545947SLisandro Dalcin         PetscInt        i, j, len, nnz, cnt, *iia = NULL, *jja = NULL;
13621a83f524SJed Brown         const PetscInt *row;
136392bb6962SLisandro Dalcin         nnz = 0;
136492bb6962SLisandro Dalcin         for (i = 0; i < na; i++) { /* count number of nonzeros */
136592bb6962SLisandro Dalcin           len = ia[i + 1] - ia[i];
136692bb6962SLisandro Dalcin           row = ja + ia[i];
136792bb6962SLisandro Dalcin           for (j = 0; j < len; j++) {
136892bb6962SLisandro Dalcin             if (row[j] == i) { /* don't count diagonal */
13699371c9d4SSatish Balay               len--;
13709371c9d4SSatish Balay               break;
137192bb6962SLisandro Dalcin             }
137292bb6962SLisandro Dalcin           }
137392bb6962SLisandro Dalcin           nnz += len;
137492bb6962SLisandro Dalcin         }
13759566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(na + 1, &iia));
13769566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(nnz, &jja));
137792bb6962SLisandro Dalcin         nnz    = 0;
137892bb6962SLisandro Dalcin         iia[0] = 0;
137992bb6962SLisandro Dalcin         for (i = 0; i < na; i++) { /* fill adjacency */
138092bb6962SLisandro Dalcin           cnt = 0;
138192bb6962SLisandro Dalcin           len = ia[i + 1] - ia[i];
138292bb6962SLisandro Dalcin           row = ja + ia[i];
138392bb6962SLisandro Dalcin           for (j = 0; j < len; j++) {
138492bb6962SLisandro Dalcin             if (row[j] != i) { /* if not diagonal */
138592bb6962SLisandro Dalcin               jja[nnz + cnt++] = row[j];
138692bb6962SLisandro Dalcin             }
138792bb6962SLisandro Dalcin           }
138892bb6962SLisandro Dalcin           nnz += cnt;
138992bb6962SLisandro Dalcin           iia[i + 1] = nnz;
139092bb6962SLisandro Dalcin         }
139192bb6962SLisandro Dalcin         /* Partitioning of the adjacency matrix */
13929566063dSJacob Faibussowitsch         PetscCall(MatCreateMPIAdj(PETSC_COMM_SELF, na, na, iia, jja, NULL, &adj));
13939566063dSJacob Faibussowitsch         PetscCall(MatPartitioningSetAdjacency(mpart, adj));
13949566063dSJacob Faibussowitsch         PetscCall(MatPartitioningSetNParts(mpart, n));
13959566063dSJacob Faibussowitsch         PetscCall(MatPartitioningApply(mpart, &ispart));
13969566063dSJacob Faibussowitsch         PetscCall(ISPartitioningToNumbering(ispart, &isnumb));
13979566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&adj));
139892bb6962SLisandro Dalcin         foundpart = PETSC_TRUE;
139992bb6962SLisandro Dalcin       }
14009566063dSJacob Faibussowitsch       PetscCall(MatRestoreRowIJ(Ad, 0, PETSC_TRUE, isbaij, &na, &ia, &ja, &done));
140192bb6962SLisandro Dalcin     }
14029566063dSJacob Faibussowitsch     PetscCall(MatPartitioningDestroy(&mpart));
140392bb6962SLisandro Dalcin   }
140492bb6962SLisandro Dalcin 
14059566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n, &is));
140692bb6962SLisandro Dalcin   *outis = is;
140792bb6962SLisandro Dalcin 
140892bb6962SLisandro Dalcin   if (!foundpart) {
140992bb6962SLisandro Dalcin     /* Partitioning by contiguous chunks of rows */
141092bb6962SLisandro Dalcin 
141192bb6962SLisandro Dalcin     PetscInt mbs   = (rend - rstart) / bs;
141292bb6962SLisandro Dalcin     PetscInt start = rstart;
141392bb6962SLisandro Dalcin     for (i = 0; i < n; i++) {
141492bb6962SLisandro Dalcin       PetscInt count = (mbs / n + ((mbs % n) > i)) * bs;
14159566063dSJacob Faibussowitsch       PetscCall(ISCreateStride(PETSC_COMM_SELF, count, start, 1, &is[i]));
141692bb6962SLisandro Dalcin       start += count;
141792bb6962SLisandro Dalcin     }
141892bb6962SLisandro Dalcin 
141992bb6962SLisandro Dalcin   } else {
142092bb6962SLisandro Dalcin     /* Partitioning by adjacency of diagonal block  */
142192bb6962SLisandro Dalcin 
142292bb6962SLisandro Dalcin     const PetscInt *numbering;
142392bb6962SLisandro Dalcin     PetscInt       *count, nidx, *indices, *newidx, start = 0;
142492bb6962SLisandro Dalcin     /* Get node count in each partition */
14259566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &count));
14269566063dSJacob Faibussowitsch     PetscCall(ISPartitioningCount(ispart, n, count));
142792bb6962SLisandro Dalcin     if (isbaij && bs > 1) { /* adjust for the block-aij case */
142892bb6962SLisandro Dalcin       for (i = 0; i < n; i++) count[i] *= bs;
142992bb6962SLisandro Dalcin     }
143092bb6962SLisandro Dalcin     /* Build indices from node numbering */
14319566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(isnumb, &nidx));
14329566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nidx, &indices));
143392bb6962SLisandro Dalcin     for (i = 0; i < nidx; i++) indices[i] = i; /* needs to be initialized */
14349566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(isnumb, &numbering));
14359566063dSJacob Faibussowitsch     PetscCall(PetscSortIntWithPermutation(nidx, numbering, indices));
14369566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(isnumb, &numbering));
143792bb6962SLisandro Dalcin     if (isbaij && bs > 1) { /* adjust for the block-aij case */
14389566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(nidx * bs, &newidx));
14392fa5cd67SKarl Rupp       for (i = 0; i < nidx; i++) {
14402fa5cd67SKarl Rupp         for (j = 0; j < bs; j++) newidx[i * bs + j] = indices[i] * bs + j;
14412fa5cd67SKarl Rupp       }
14429566063dSJacob Faibussowitsch       PetscCall(PetscFree(indices));
144392bb6962SLisandro Dalcin       nidx *= bs;
144492bb6962SLisandro Dalcin       indices = newidx;
144592bb6962SLisandro Dalcin     }
144692bb6962SLisandro Dalcin     /* Shift to get global indices */
144792bb6962SLisandro Dalcin     for (i = 0; i < nidx; i++) indices[i] += rstart;
144892bb6962SLisandro Dalcin 
144992bb6962SLisandro Dalcin     /* Build the index sets for each block */
145092bb6962SLisandro Dalcin     for (i = 0; i < n; i++) {
14519566063dSJacob Faibussowitsch       PetscCall(ISCreateGeneral(PETSC_COMM_SELF, count[i], &indices[start], PETSC_COPY_VALUES, &is[i]));
14529566063dSJacob Faibussowitsch       PetscCall(ISSort(is[i]));
145392bb6962SLisandro Dalcin       start += count[i];
145492bb6962SLisandro Dalcin     }
145592bb6962SLisandro Dalcin 
14569566063dSJacob Faibussowitsch     PetscCall(PetscFree(count));
14579566063dSJacob Faibussowitsch     PetscCall(PetscFree(indices));
14589566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&isnumb));
14599566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&ispart));
146092bb6962SLisandro Dalcin   }
146192bb6962SLisandro Dalcin   PetscFunctionReturn(0);
146292bb6962SLisandro Dalcin }
146392bb6962SLisandro Dalcin 
146492bb6962SLisandro Dalcin /*@C
146592bb6962SLisandro Dalcin    PCASMDestroySubdomains - Destroys the index sets created with
146692bb6962SLisandro Dalcin    PCASMCreateSubdomains(). Should be called after setting subdomains
146792bb6962SLisandro Dalcin    with PCASMSetLocalSubdomains().
146892bb6962SLisandro Dalcin 
146992bb6962SLisandro Dalcin    Collective
147092bb6962SLisandro Dalcin 
147192bb6962SLisandro Dalcin    Input Parameters:
147292bb6962SLisandro Dalcin +  n - the number of index sets
14732b691e39SMatthew Knepley .  is - the array of index sets
14740298fd71SBarry Smith -  is_local - the array of local index sets, can be NULL
147592bb6962SLisandro Dalcin 
147692bb6962SLisandro Dalcin    Level: advanced
147792bb6962SLisandro Dalcin 
1478db781477SPatrick Sanan .seealso: `PCASMCreateSubdomains()`, `PCASMSetLocalSubdomains()`
147992bb6962SLisandro Dalcin @*/
14809371c9d4SSatish Balay PetscErrorCode PCASMDestroySubdomains(PetscInt n, IS is[], IS is_local[]) {
148192bb6962SLisandro Dalcin   PetscInt i;
14825fd66863SKarl Rupp 
148392bb6962SLisandro Dalcin   PetscFunctionBegin;
1484a35b7b57SDmitry Karpeev   if (n <= 0) PetscFunctionReturn(0);
1485a35b7b57SDmitry Karpeev   if (is) {
148692bb6962SLisandro Dalcin     PetscValidPointer(is, 2);
14879566063dSJacob Faibussowitsch     for (i = 0; i < n; i++) PetscCall(ISDestroy(&is[i]));
14889566063dSJacob Faibussowitsch     PetscCall(PetscFree(is));
1489a35b7b57SDmitry Karpeev   }
14902b691e39SMatthew Knepley   if (is_local) {
14912b691e39SMatthew Knepley     PetscValidPointer(is_local, 3);
14929566063dSJacob Faibussowitsch     for (i = 0; i < n; i++) PetscCall(ISDestroy(&is_local[i]));
14939566063dSJacob Faibussowitsch     PetscCall(PetscFree(is_local));
14942b691e39SMatthew Knepley   }
149592bb6962SLisandro Dalcin   PetscFunctionReturn(0);
149692bb6962SLisandro Dalcin }
149792bb6962SLisandro Dalcin 
14984b9ad928SBarry Smith /*@
14994b9ad928SBarry Smith    PCASMCreateSubdomains2D - Creates the index sets for the overlapping Schwarz
15004b9ad928SBarry Smith    preconditioner for a two-dimensional problem on a regular grid.
15014b9ad928SBarry Smith 
15024b9ad928SBarry Smith    Not Collective
15034b9ad928SBarry Smith 
15044b9ad928SBarry Smith    Input Parameters:
15056b867d5aSJose E. Roman +  m   - the number of mesh points in the x direction
15066b867d5aSJose E. Roman .  n   - the number of mesh points in the y direction
15076b867d5aSJose E. Roman .  M   - the number of subdomains in the x direction
15086b867d5aSJose E. Roman .  N   - the number of subdomains in the y direction
15094b9ad928SBarry Smith .  dof - degrees of freedom per node
15104b9ad928SBarry Smith -  overlap - overlap in mesh lines
15114b9ad928SBarry Smith 
15124b9ad928SBarry Smith    Output Parameters:
15134b9ad928SBarry Smith +  Nsub - the number of subdomains created
15143d03bb48SJed Brown .  is - array of index sets defining overlapping (if overlap > 0) subdomains
15153d03bb48SJed Brown -  is_local - array of index sets defining non-overlapping subdomains
15164b9ad928SBarry Smith 
15174b9ad928SBarry Smith    Note:
15184b9ad928SBarry Smith    Presently PCAMSCreateSubdomains2d() is valid only for sequential
15194b9ad928SBarry Smith    preconditioners.  More general related routines are
15204b9ad928SBarry Smith    PCASMSetTotalSubdomains() and PCASMSetLocalSubdomains().
15214b9ad928SBarry Smith 
15224b9ad928SBarry Smith    Level: advanced
15234b9ad928SBarry Smith 
1524db781477SPatrick Sanan .seealso: `PCASMSetTotalSubdomains()`, `PCASMSetLocalSubdomains()`, `PCASMGetSubKSP()`,
1525db781477SPatrick Sanan           `PCASMSetOverlap()`
15264b9ad928SBarry Smith @*/
15279371c9d4SSatish Balay PetscErrorCode PCASMCreateSubdomains2D(PetscInt m, PetscInt n, PetscInt M, PetscInt N, PetscInt dof, PetscInt overlap, PetscInt *Nsub, IS **is, IS **is_local) {
15283d03bb48SJed Brown   PetscInt i, j, height, width, ystart, xstart, yleft, yright, xleft, xright, loc_outer;
152913f74950SBarry Smith   PetscInt nidx, *idx, loc, ii, jj, count;
15304b9ad928SBarry Smith 
15314b9ad928SBarry Smith   PetscFunctionBegin;
15327827d75bSBarry Smith   PetscCheck(dof == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "dof must be 1");
15334b9ad928SBarry Smith 
15344b9ad928SBarry Smith   *Nsub = N * M;
15359566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(*Nsub, is));
15369566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(*Nsub, is_local));
15374b9ad928SBarry Smith   ystart    = 0;
15383d03bb48SJed Brown   loc_outer = 0;
15394b9ad928SBarry Smith   for (i = 0; i < N; i++) {
15404b9ad928SBarry Smith     height = n / N + ((n % N) > i); /* height of subdomain */
15417827d75bSBarry Smith     PetscCheck(height >= 2, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many N subdomains for mesh dimension n");
15429371c9d4SSatish Balay     yleft = ystart - overlap;
15439371c9d4SSatish Balay     if (yleft < 0) yleft = 0;
15449371c9d4SSatish Balay     yright = ystart + height + overlap;
15459371c9d4SSatish Balay     if (yright > n) yright = n;
15464b9ad928SBarry Smith     xstart = 0;
15474b9ad928SBarry Smith     for (j = 0; j < M; j++) {
15484b9ad928SBarry Smith       width = m / M + ((m % M) > j); /* width of subdomain */
15497827d75bSBarry Smith       PetscCheck(width >= 2, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many M subdomains for mesh dimension m");
15509371c9d4SSatish Balay       xleft = xstart - overlap;
15519371c9d4SSatish Balay       if (xleft < 0) xleft = 0;
15529371c9d4SSatish Balay       xright = xstart + width + overlap;
15539371c9d4SSatish Balay       if (xright > m) xright = m;
15544b9ad928SBarry Smith       nidx = (xright - xleft) * (yright - yleft);
15559566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(nidx, &idx));
15564b9ad928SBarry Smith       loc = 0;
15574b9ad928SBarry Smith       for (ii = yleft; ii < yright; ii++) {
15584b9ad928SBarry Smith         count = m * ii + xleft;
15592fa5cd67SKarl Rupp         for (jj = xleft; jj < xright; jj++) idx[loc++] = count++;
15604b9ad928SBarry Smith       }
15619566063dSJacob Faibussowitsch       PetscCall(ISCreateGeneral(PETSC_COMM_SELF, nidx, idx, PETSC_COPY_VALUES, (*is) + loc_outer));
15623d03bb48SJed Brown       if (overlap == 0) {
15639566063dSJacob Faibussowitsch         PetscCall(PetscObjectReference((PetscObject)(*is)[loc_outer]));
15642fa5cd67SKarl Rupp 
15653d03bb48SJed Brown         (*is_local)[loc_outer] = (*is)[loc_outer];
15663d03bb48SJed Brown       } else {
15673d03bb48SJed Brown         for (loc = 0, ii = ystart; ii < ystart + height; ii++) {
15689371c9d4SSatish Balay           for (jj = xstart; jj < xstart + width; jj++) { idx[loc++] = m * ii + jj; }
15693d03bb48SJed Brown         }
15709566063dSJacob Faibussowitsch         PetscCall(ISCreateGeneral(PETSC_COMM_SELF, loc, idx, PETSC_COPY_VALUES, *is_local + loc_outer));
15713d03bb48SJed Brown       }
15729566063dSJacob Faibussowitsch       PetscCall(PetscFree(idx));
15734b9ad928SBarry Smith       xstart += width;
15743d03bb48SJed Brown       loc_outer++;
15754b9ad928SBarry Smith     }
15764b9ad928SBarry Smith     ystart += height;
15774b9ad928SBarry Smith   }
15789566063dSJacob Faibussowitsch   for (i = 0; i < *Nsub; i++) PetscCall(ISSort((*is)[i]));
15794b9ad928SBarry Smith   PetscFunctionReturn(0);
15804b9ad928SBarry Smith }
15814b9ad928SBarry Smith 
15824b9ad928SBarry Smith /*@C
15834b9ad928SBarry Smith     PCASMGetLocalSubdomains - Gets the local subdomains (for this processor
15844b9ad928SBarry Smith     only) for the additive Schwarz preconditioner.
15854b9ad928SBarry Smith 
1586ad4df100SBarry Smith     Not Collective
15874b9ad928SBarry Smith 
15884b9ad928SBarry Smith     Input Parameter:
15894b9ad928SBarry Smith .   pc - the preconditioner context
15904b9ad928SBarry Smith 
15914b9ad928SBarry Smith     Output Parameters:
1592f17a6784SPierre Jolivet +   n - if requested, the number of subdomains for this processor (default value = 1)
1593f17a6784SPierre Jolivet .   is - if requested, the index sets that define the subdomains for this processor
1594f17a6784SPierre Jolivet -   is_local - if requested, the index sets that define the local part of the subdomains for this processor (can be NULL)
15954b9ad928SBarry Smith 
15964b9ad928SBarry Smith     Notes:
15974b9ad928SBarry Smith     The IS numbering is in the parallel, global numbering of the vector.
15984b9ad928SBarry Smith 
15994b9ad928SBarry Smith     Level: advanced
16004b9ad928SBarry Smith 
1601db781477SPatrick Sanan .seealso: `PCASMSetTotalSubdomains()`, `PCASMSetOverlap()`, `PCASMGetSubKSP()`,
1602db781477SPatrick Sanan           `PCASMCreateSubdomains2D()`, `PCASMSetLocalSubdomains()`, `PCASMGetLocalSubmatrices()`
16034b9ad928SBarry Smith @*/
16049371c9d4SSatish Balay PetscErrorCode PCASMGetLocalSubdomains(PC pc, PetscInt *n, IS *is[], IS *is_local[]) {
16052a808120SBarry Smith   PC_ASM   *osm = (PC_ASM *)pc->data;
1606ace3abfcSBarry Smith   PetscBool match;
16074b9ad928SBarry Smith 
16084b9ad928SBarry Smith   PetscFunctionBegin;
16090700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1610f17a6784SPierre Jolivet   if (n) PetscValidIntPointer(n, 2);
161192bb6962SLisandro Dalcin   if (is) PetscValidPointer(is, 3);
1612f17a6784SPierre Jolivet   if (is_local) PetscValidPointer(is_local, 4);
16139566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCASM, &match));
161428b400f6SJacob Faibussowitsch   PetscCheck(match, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "PC is not a PCASM");
16154b9ad928SBarry Smith   if (n) *n = osm->n_local_true;
16164b9ad928SBarry Smith   if (is) *is = osm->is;
16172b691e39SMatthew Knepley   if (is_local) *is_local = osm->is_local;
16184b9ad928SBarry Smith   PetscFunctionReturn(0);
16194b9ad928SBarry Smith }
16204b9ad928SBarry Smith 
16214b9ad928SBarry Smith /*@C
16224b9ad928SBarry Smith     PCASMGetLocalSubmatrices - Gets the local submatrices (for this processor
16234b9ad928SBarry Smith     only) for the additive Schwarz preconditioner.
16244b9ad928SBarry Smith 
1625ad4df100SBarry Smith     Not Collective
16264b9ad928SBarry Smith 
16274b9ad928SBarry Smith     Input Parameter:
16284b9ad928SBarry Smith .   pc - the preconditioner context
16294b9ad928SBarry Smith 
16304b9ad928SBarry Smith     Output Parameters:
1631f17a6784SPierre Jolivet +   n - if requested, the number of matrices for this processor (default value = 1)
1632f17a6784SPierre Jolivet -   mat - if requested, the matrices
16334b9ad928SBarry Smith 
16344b9ad928SBarry Smith     Level: advanced
16354b9ad928SBarry Smith 
163695452b02SPatrick Sanan     Notes:
163734a84908Sprj-     Call after PCSetUp() (or KSPSetUp()) but before PCApply() and before PCSetUpOnBlocks())
1638cf739d55SBarry Smith 
163934a84908Sprj-            Usually one would use PCSetModifySubMatrices() to change the submatrices in building the preconditioner.
1640cf739d55SBarry Smith 
1641db781477SPatrick Sanan .seealso: `PCASMSetTotalSubdomains()`, `PCASMSetOverlap()`, `PCASMGetSubKSP()`,
1642db781477SPatrick Sanan           `PCASMCreateSubdomains2D()`, `PCASMSetLocalSubdomains()`, `PCASMGetLocalSubdomains()`, `PCSetModifySubMatrices()`
16434b9ad928SBarry Smith @*/
16449371c9d4SSatish Balay PetscErrorCode PCASMGetLocalSubmatrices(PC pc, PetscInt *n, Mat *mat[]) {
16454b9ad928SBarry Smith   PC_ASM   *osm;
1646ace3abfcSBarry Smith   PetscBool match;
16474b9ad928SBarry Smith 
16484b9ad928SBarry Smith   PetscFunctionBegin;
16490700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1650f17a6784SPierre Jolivet   if (n) PetscValidIntPointer(n, 2);
165192bb6962SLisandro Dalcin   if (mat) PetscValidPointer(mat, 3);
165228b400f6SJacob Faibussowitsch   PetscCheck(pc->setupcalled, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must call after KSPSetUp() or PCSetUp().");
16539566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCASM, &match));
165492bb6962SLisandro Dalcin   if (!match) {
165592bb6962SLisandro Dalcin     if (n) *n = 0;
16560298fd71SBarry Smith     if (mat) *mat = NULL;
165792bb6962SLisandro Dalcin   } else {
16584b9ad928SBarry Smith     osm = (PC_ASM *)pc->data;
16594b9ad928SBarry Smith     if (n) *n = osm->n_local_true;
16604b9ad928SBarry Smith     if (mat) *mat = osm->pmat;
166192bb6962SLisandro Dalcin   }
16624b9ad928SBarry Smith   PetscFunctionReturn(0);
16634b9ad928SBarry Smith }
1664d709ea83SDmitry Karpeev 
1665d709ea83SDmitry Karpeev /*@
1666d709ea83SDmitry Karpeev     PCASMSetDMSubdomains - Indicates whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible.
1667f1ee410cSBarry Smith 
1668d709ea83SDmitry Karpeev     Logically Collective
1669d709ea83SDmitry Karpeev 
1670d8d19677SJose E. Roman     Input Parameters:
1671d709ea83SDmitry Karpeev +   pc  - the preconditioner
1672d709ea83SDmitry Karpeev -   flg - boolean indicating whether to use subdomains defined by the DM
1673d709ea83SDmitry Karpeev 
1674d709ea83SDmitry Karpeev     Options Database Key:
1675147403d9SBarry Smith .   -pc_asm_dm_subdomains <bool> - use subdomains defined by the DM
1676d709ea83SDmitry Karpeev 
1677d709ea83SDmitry Karpeev     Level: intermediate
1678d709ea83SDmitry Karpeev 
1679d709ea83SDmitry Karpeev     Notes:
1680d709ea83SDmitry Karpeev     PCASMSetTotalSubdomains() and PCASMSetOverlap() take precedence over PCASMSetDMSubdomains(),
1681d709ea83SDmitry Karpeev     so setting either of the first two effectively turns the latter off.
1682d709ea83SDmitry Karpeev 
1683db781477SPatrick Sanan .seealso: `PCASMGetDMSubdomains()`, `PCASMSetTotalSubdomains()`, `PCASMSetOverlap()`
1684db781477SPatrick Sanan           `PCASMCreateSubdomains2D()`, `PCASMSetLocalSubdomains()`, `PCASMGetLocalSubdomains()`
1685d709ea83SDmitry Karpeev @*/
16869371c9d4SSatish Balay PetscErrorCode PCASMSetDMSubdomains(PC pc, PetscBool flg) {
1687d709ea83SDmitry Karpeev   PC_ASM   *osm = (PC_ASM *)pc->data;
1688d709ea83SDmitry Karpeev   PetscBool match;
1689d709ea83SDmitry Karpeev 
1690d709ea83SDmitry Karpeev   PetscFunctionBegin;
1691d709ea83SDmitry Karpeev   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1692d709ea83SDmitry Karpeev   PetscValidLogicalCollectiveBool(pc, flg, 2);
169328b400f6SJacob Faibussowitsch   PetscCheck(!pc->setupcalled, ((PetscObject)pc)->comm, PETSC_ERR_ARG_WRONGSTATE, "Not for a setup PC.");
16949566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCASM, &match));
16959371c9d4SSatish Balay   if (match) { osm->dm_subdomains = flg; }
1696d709ea83SDmitry Karpeev   PetscFunctionReturn(0);
1697d709ea83SDmitry Karpeev }
1698d709ea83SDmitry Karpeev 
1699d709ea83SDmitry Karpeev /*@
1700d709ea83SDmitry Karpeev     PCASMGetDMSubdomains - Returns flag indicating whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible.
1701d709ea83SDmitry Karpeev     Not Collective
1702d709ea83SDmitry Karpeev 
1703d709ea83SDmitry Karpeev     Input Parameter:
1704d709ea83SDmitry Karpeev .   pc  - the preconditioner
1705d709ea83SDmitry Karpeev 
1706d709ea83SDmitry Karpeev     Output Parameter:
1707d709ea83SDmitry Karpeev .   flg - boolean indicating whether to use subdomains defined by the DM
1708d709ea83SDmitry Karpeev 
1709d709ea83SDmitry Karpeev     Level: intermediate
1710d709ea83SDmitry Karpeev 
1711db781477SPatrick Sanan .seealso: `PCASMSetDMSubdomains()`, `PCASMSetTotalSubdomains()`, `PCASMSetOverlap()`
1712db781477SPatrick Sanan           `PCASMCreateSubdomains2D()`, `PCASMSetLocalSubdomains()`, `PCASMGetLocalSubdomains()`
1713d709ea83SDmitry Karpeev @*/
17149371c9d4SSatish Balay PetscErrorCode PCASMGetDMSubdomains(PC pc, PetscBool *flg) {
1715d709ea83SDmitry Karpeev   PC_ASM   *osm = (PC_ASM *)pc->data;
1716d709ea83SDmitry Karpeev   PetscBool match;
1717d709ea83SDmitry Karpeev 
1718d709ea83SDmitry Karpeev   PetscFunctionBegin;
1719d709ea83SDmitry Karpeev   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1720534a8f05SLisandro Dalcin   PetscValidBoolPointer(flg, 2);
17219566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCASM, &match));
172256ea09ceSStefano Zampini   if (match) *flg = osm->dm_subdomains;
172356ea09ceSStefano Zampini   else *flg = PETSC_FALSE;
1724d709ea83SDmitry Karpeev   PetscFunctionReturn(0);
1725d709ea83SDmitry Karpeev }
172680ec0b7dSPatrick Sanan 
172780ec0b7dSPatrick Sanan /*@
172880ec0b7dSPatrick Sanan      PCASMGetSubMatType - Gets the matrix type used for ASM subsolves, as a string.
172980ec0b7dSPatrick Sanan 
173080ec0b7dSPatrick Sanan    Not Collective
173180ec0b7dSPatrick Sanan 
173280ec0b7dSPatrick Sanan    Input Parameter:
173380ec0b7dSPatrick Sanan .  pc - the PC
173480ec0b7dSPatrick Sanan 
173580ec0b7dSPatrick Sanan    Output Parameter:
1736f1ee410cSBarry Smith .  -pc_asm_sub_mat_type - name of matrix type
173780ec0b7dSPatrick Sanan 
173880ec0b7dSPatrick Sanan    Level: advanced
173980ec0b7dSPatrick Sanan 
1740db781477SPatrick Sanan .seealso: `PCASMSetSubMatType()`, `PCASM`, `PCSetType()`, `VecSetType()`, `MatType`, `Mat`
174180ec0b7dSPatrick Sanan @*/
17429371c9d4SSatish Balay PetscErrorCode PCASMGetSubMatType(PC pc, MatType *sub_mat_type) {
174356ea09ceSStefano Zampini   PetscFunctionBegin;
174456ea09ceSStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1745cac4c232SBarry Smith   PetscTryMethod(pc, "PCASMGetSubMatType_C", (PC, MatType *), (pc, sub_mat_type));
174680ec0b7dSPatrick Sanan   PetscFunctionReturn(0);
174780ec0b7dSPatrick Sanan }
174880ec0b7dSPatrick Sanan 
174980ec0b7dSPatrick Sanan /*@
175080ec0b7dSPatrick Sanan      PCASMSetSubMatType - Set the type of matrix used for ASM subsolves
175180ec0b7dSPatrick Sanan 
175280ec0b7dSPatrick Sanan    Collective on Mat
175380ec0b7dSPatrick Sanan 
175480ec0b7dSPatrick Sanan    Input Parameters:
175580ec0b7dSPatrick Sanan +  pc             - the PC object
175680ec0b7dSPatrick Sanan -  sub_mat_type   - matrix type
175780ec0b7dSPatrick Sanan 
175880ec0b7dSPatrick Sanan    Options Database Key:
175980ec0b7dSPatrick 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.
176080ec0b7dSPatrick Sanan 
176180ec0b7dSPatrick Sanan    Notes:
176280ec0b7dSPatrick Sanan    See "${PETSC_DIR}/include/petscmat.h" for available types
176380ec0b7dSPatrick Sanan 
176480ec0b7dSPatrick Sanan   Level: advanced
176580ec0b7dSPatrick Sanan 
1766db781477SPatrick Sanan .seealso: `PCASMGetSubMatType()`, `PCASM`, `PCSetType()`, `VecSetType()`, `MatType`, `Mat`
176780ec0b7dSPatrick Sanan @*/
17689371c9d4SSatish Balay PetscErrorCode PCASMSetSubMatType(PC pc, MatType sub_mat_type) {
176956ea09ceSStefano Zampini   PetscFunctionBegin;
177056ea09ceSStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1771cac4c232SBarry Smith   PetscTryMethod(pc, "PCASMSetSubMatType_C", (PC, MatType), (pc, sub_mat_type));
177280ec0b7dSPatrick Sanan   PetscFunctionReturn(0);
177380ec0b7dSPatrick Sanan }
1774