xref: /petsc/src/ksp/pc/impls/asm/asm.c (revision 4dbf25a8fa98e38799e7b47dcb2d8a9309975f41)
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*/
1446233b44SBarry Smith #include <petsc/private/matimpl.h>
154b9ad928SBarry Smith 
16d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCView_ASM(PC pc, PetscViewer viewer)
17d71ae5a4SJacob Faibussowitsch {
1892bb6962SLisandro Dalcin   PC_ASM           *osm = (PC_ASM *)pc->data;
1913f74950SBarry Smith   PetscMPIInt       rank;
204d219a6aSLisandro Dalcin   PetscInt          i, bsz;
21ace3abfcSBarry Smith   PetscBool         iascii, isstring;
224b9ad928SBarry Smith   PetscViewer       sviewer;
23ed155784SPierre Jolivet   PetscViewerFormat format;
24ed155784SPierre Jolivet   const char       *prefix;
254b9ad928SBarry Smith 
264b9ad928SBarry Smith   PetscFunctionBegin;
279566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
289566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSTRING, &isstring));
2932077d6dSBarry Smith   if (iascii) {
303d03bb48SJed Brown     char overlaps[256] = "user-defined overlap", blocks[256] = "total subdomain blocks not yet set";
3163a3b9bcSJacob Faibussowitsch     if (osm->overlap >= 0) PetscCall(PetscSNPrintf(overlaps, sizeof(overlaps), "amount of overlap = %" PetscInt_FMT, osm->overlap));
3263a3b9bcSJacob Faibussowitsch     if (osm->n > 0) PetscCall(PetscSNPrintf(blocks, sizeof(blocks), "total subdomain blocks = %" PetscInt_FMT, osm->n));
339566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  %s, %s\n", blocks, overlaps));
349566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "  restriction/interpolation type - %s\n", PCASMTypes[osm->type]));
359566063dSJacob Faibussowitsch     if (osm->dm_subdomains) PetscCall(PetscViewerASCIIPrintf(viewer, "  Additive Schwarz: using DM to define subdomains\n"));
369566063dSJacob Faibussowitsch     if (osm->loctype != PC_COMPOSITE_ADDITIVE) PetscCall(PetscViewerASCIIPrintf(viewer, "  Additive Schwarz: local solve composition type - %s\n", PCCompositeTypes[osm->loctype]));
379566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank));
389566063dSJacob Faibussowitsch     PetscCall(PetscViewerGetFormat(viewer, &format));
39ed155784SPierre Jolivet     if (format != PETSC_VIEWER_ASCII_INFO_DETAIL) {
404d219a6aSLisandro Dalcin       if (osm->ksp) {
419566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  Local solver information for first block is in the following KSP and PC objects on rank 0:\n"));
429566063dSJacob Faibussowitsch         PetscCall(PCGetOptionsPrefix(pc, &prefix));
439566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPrintf(viewer, "  Use -%sksp_view ::ascii_info_detail to display information for all blocks\n", prefix ? prefix : ""));
449566063dSJacob Faibussowitsch         PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
45dd400576SPatrick Sanan         if (rank == 0) {
46b4025f61SBarry Smith           PetscCall(PetscViewerASCIIPushTab(sviewer));
479566063dSJacob Faibussowitsch           PetscCall(KSPView(osm->ksp[0], sviewer));
48b4025f61SBarry Smith           PetscCall(PetscViewerASCIIPopTab(sviewer));
494b9ad928SBarry Smith         }
509566063dSJacob Faibussowitsch         PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
514d219a6aSLisandro Dalcin       }
524b9ad928SBarry Smith     } else {
539566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushSynchronized(viewer));
54835f2295SStefano Zampini       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "  [%d] number of local blocks = %" PetscInt_FMT "\n", rank, osm->n_local_true));
559566063dSJacob Faibussowitsch       PetscCall(PetscViewerFlush(viewer));
569566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "  Local solver information for each block is in the following KSP and PC objects:\n"));
579566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushTab(viewer));
589566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPrintf(viewer, "- - - - - - - - - - - - - - - - - -\n"));
599566063dSJacob Faibussowitsch       PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
60dd2fa690SBarry Smith       for (i = 0; i < osm->n_local_true; i++) {
619566063dSJacob Faibussowitsch         PetscCall(ISGetLocalSize(osm->is[i], &bsz));
62835f2295SStefano Zampini         PetscCall(PetscViewerASCIIPrintf(sviewer, "[%d] local block number %" PetscInt_FMT ", size = %" PetscInt_FMT "\n", rank, i, bsz));
639566063dSJacob Faibussowitsch         PetscCall(KSPView(osm->ksp[i], sviewer));
64b4025f61SBarry Smith         PetscCall(PetscViewerASCIIPrintf(sviewer, "- - - - - - - - - - - - - - - - - -\n"));
654b9ad928SBarry Smith       }
669566063dSJacob Faibussowitsch       PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer));
679566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopTab(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   }
763ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
774b9ad928SBarry Smith }
784b9ad928SBarry Smith 
79d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCASMPrintSubdomains(PC pc)
80d71ae5a4SJacob Faibussowitsch {
8192bb6962SLisandro Dalcin   PC_ASM         *osm = (PC_ASM *)pc->data;
8292bb6962SLisandro Dalcin   const char     *prefix;
8392bb6962SLisandro Dalcin   char            fname[PETSC_MAX_PATH_LEN + 1];
84643535aeSDmitry Karpeev   PetscViewer     viewer, sviewer;
85643535aeSDmitry Karpeev   char           *s;
8692bb6962SLisandro Dalcin   PetscInt        i, j, nidx;
8792bb6962SLisandro Dalcin   const PetscInt *idx;
88643535aeSDmitry Karpeev   PetscMPIInt     rank, size;
89846783a0SBarry Smith 
9092bb6962SLisandro Dalcin   PetscFunctionBegin;
919566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size));
929566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank));
939566063dSJacob Faibussowitsch   PetscCall(PCGetOptionsPrefix(pc, &prefix));
949566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetString(NULL, prefix, "-pc_asm_print_subdomains", fname, sizeof(fname), NULL));
95c6a7a370SJeremy L Thompson   if (fname[0] == 0) PetscCall(PetscStrncpy(fname, "stdout", sizeof(fname)));
969566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIOpen(PetscObjectComm((PetscObject)pc), fname, &viewer));
97643535aeSDmitry Karpeev   for (i = 0; i < osm->n_local; i++) {
98643535aeSDmitry Karpeev     if (i < osm->n_local_true) {
999566063dSJacob Faibussowitsch       PetscCall(ISGetLocalSize(osm->is[i], &nidx));
1009566063dSJacob Faibussowitsch       PetscCall(ISGetIndices(osm->is[i], &idx));
101643535aeSDmitry Karpeev       /* Print to a string viewer; no more than 15 characters per index plus 512 char for the header.*/
10236a9e3b9SBarry Smith #define len 16 * (nidx + 1) + 512
1039566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(len, &s));
1049566063dSJacob Faibussowitsch       PetscCall(PetscViewerStringOpen(PETSC_COMM_SELF, s, len, &sviewer));
10536a9e3b9SBarry Smith #undef len
10663a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerStringSPrintf(sviewer, "[%d:%d] Subdomain %" PetscInt_FMT " with overlap:\n", rank, size, i));
10748a46eb9SPierre Jolivet       for (j = 0; j < nidx; j++) PetscCall(PetscViewerStringSPrintf(sviewer, "%" PetscInt_FMT " ", idx[j]));
1089566063dSJacob Faibussowitsch       PetscCall(ISRestoreIndices(osm->is[i], &idx));
1099566063dSJacob Faibussowitsch       PetscCall(PetscViewerStringSPrintf(sviewer, "\n"));
1109566063dSJacob Faibussowitsch       PetscCall(PetscViewerDestroy(&sviewer));
1119566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushSynchronized(viewer));
11263a3b9bcSJacob Faibussowitsch       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "%s", s));
1139566063dSJacob Faibussowitsch       PetscCall(PetscViewerFlush(viewer));
1149566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopSynchronized(viewer));
1159566063dSJacob Faibussowitsch       PetscCall(PetscFree(s));
1162b691e39SMatthew Knepley       if (osm->is_local) {
117643535aeSDmitry Karpeev         /* Print to a string viewer; no more than 15 characters per index plus 512 char for the header.*/
11836a9e3b9SBarry Smith #define len 16 * (nidx + 1) + 512
1199566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(len, &s));
1209566063dSJacob Faibussowitsch         PetscCall(PetscViewerStringOpen(PETSC_COMM_SELF, s, len, &sviewer));
12136a9e3b9SBarry Smith #undef len
12263a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerStringSPrintf(sviewer, "[%d:%d] Subdomain %" PetscInt_FMT " without overlap:\n", rank, size, i));
1239566063dSJacob Faibussowitsch         PetscCall(ISGetLocalSize(osm->is_local[i], &nidx));
1249566063dSJacob Faibussowitsch         PetscCall(ISGetIndices(osm->is_local[i], &idx));
12548a46eb9SPierre Jolivet         for (j = 0; j < nidx; j++) PetscCall(PetscViewerStringSPrintf(sviewer, "%" PetscInt_FMT " ", idx[j]));
1269566063dSJacob Faibussowitsch         PetscCall(ISRestoreIndices(osm->is_local[i], &idx));
1279566063dSJacob Faibussowitsch         PetscCall(PetscViewerStringSPrintf(sviewer, "\n"));
1289566063dSJacob Faibussowitsch         PetscCall(PetscViewerDestroy(&sviewer));
1299566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPushSynchronized(viewer));
13063a3b9bcSJacob Faibussowitsch         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "%s", s));
1319566063dSJacob Faibussowitsch         PetscCall(PetscViewerFlush(viewer));
1329566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPopSynchronized(viewer));
1339566063dSJacob Faibussowitsch         PetscCall(PetscFree(s));
134643535aeSDmitry Karpeev       }
1352fa5cd67SKarl Rupp     } else {
136643535aeSDmitry Karpeev       /* Participate in collective viewer calls. */
1379566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPushSynchronized(viewer));
1389566063dSJacob Faibussowitsch       PetscCall(PetscViewerFlush(viewer));
1399566063dSJacob Faibussowitsch       PetscCall(PetscViewerASCIIPopSynchronized(viewer));
140643535aeSDmitry Karpeev       /* Assume either all ranks have is_local or none do. */
141643535aeSDmitry Karpeev       if (osm->is_local) {
1429566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPushSynchronized(viewer));
1439566063dSJacob Faibussowitsch         PetscCall(PetscViewerFlush(viewer));
1449566063dSJacob Faibussowitsch         PetscCall(PetscViewerASCIIPopSynchronized(viewer));
145643535aeSDmitry Karpeev       }
146fdc48646SDmitry Karpeev     }
14792bb6962SLisandro Dalcin   }
1489566063dSJacob Faibussowitsch   PetscCall(PetscViewerFlush(viewer));
1499566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&viewer));
1503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
15192bb6962SLisandro Dalcin }
15292bb6962SLisandro Dalcin 
153d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_ASM(PC pc)
154d71ae5a4SJacob Faibussowitsch {
1554b9ad928SBarry Smith   PC_ASM       *osm = (PC_ASM *)pc->data;
15657501b6eSBarry Smith   PetscBool     flg;
15787e86712SBarry Smith   PetscInt      i, m, m_local;
1584b9ad928SBarry Smith   MatReuse      scall = MAT_REUSE_MATRIX;
1594b9ad928SBarry Smith   IS            isl;
1604b9ad928SBarry Smith   KSP           ksp;
1614b9ad928SBarry Smith   PC            subpc;
1622dcb1b2aSMatthew Knepley   const char   *prefix, *pprefix;
16323ce1328SBarry Smith   Vec           vec;
1640298fd71SBarry Smith   DM           *domain_dm = NULL;
1654849c82aSBarry Smith   MatNullSpace *nullsp    = NULL;
1664b9ad928SBarry Smith 
1674b9ad928SBarry Smith   PetscFunctionBegin;
1684b9ad928SBarry Smith   if (!pc->setupcalled) {
169265a2120SBarry Smith     PetscInt m;
17092bb6962SLisandro Dalcin 
1712b837212SDmitry Karpeev     /* Note: if subdomains have been set either via PCASMSetTotalSubdomains() or via PCASMSetLocalSubdomains(), osm->n_local_true will not be PETSC_DECIDE */
1722b837212SDmitry Karpeev     if (osm->n_local_true == PETSC_DECIDE) {
17369ca1f37SDmitry Karpeev       /* no subdomains given */
17465db9045SDmitry Karpeev       /* try pc->dm first, if allowed */
175d709ea83SDmitry Karpeev       if (osm->dm_subdomains && pc->dm) {
176feb221f8SDmitry Karpeev         PetscInt num_domains, d;
177feb221f8SDmitry Karpeev         char   **domain_names;
1788d4ac253SDmitry Karpeev         IS      *inner_domain_is, *outer_domain_is;
1799566063dSJacob Faibussowitsch         PetscCall(DMCreateDomainDecomposition(pc->dm, &num_domains, &domain_names, &inner_domain_is, &outer_domain_is, &domain_dm));
180704f0395SPatrick Sanan         osm->overlap = -1; /* We do not want to increase the overlap of the IS.
181704f0395SPatrick Sanan                               A future improvement of this code might allow one to use
182704f0395SPatrick Sanan                               DM-defined subdomains and also increase the overlap,
183704f0395SPatrick Sanan                               but that is not currently supported */
1841baa6e33SBarry Smith         if (num_domains) PetscCall(PCASMSetLocalSubdomains(pc, num_domains, outer_domain_is, inner_domain_is));
185feb221f8SDmitry Karpeev         for (d = 0; d < num_domains; ++d) {
1869566063dSJacob Faibussowitsch           if (domain_names) PetscCall(PetscFree(domain_names[d]));
1879566063dSJacob Faibussowitsch           if (inner_domain_is) PetscCall(ISDestroy(&inner_domain_is[d]));
1889566063dSJacob Faibussowitsch           if (outer_domain_is) PetscCall(ISDestroy(&outer_domain_is[d]));
189feb221f8SDmitry Karpeev         }
1909566063dSJacob Faibussowitsch         PetscCall(PetscFree(domain_names));
1919566063dSJacob Faibussowitsch         PetscCall(PetscFree(inner_domain_is));
1929566063dSJacob Faibussowitsch         PetscCall(PetscFree(outer_domain_is));
193feb221f8SDmitry Karpeev       }
1942b837212SDmitry Karpeev       if (osm->n_local_true == PETSC_DECIDE) {
19569ca1f37SDmitry Karpeev         /* still no subdomains; use one subdomain per processor */
1962b837212SDmitry Karpeev         osm->n_local_true = 1;
197feb221f8SDmitry Karpeev       }
1982b837212SDmitry Karpeev     }
1992b837212SDmitry Karpeev     { /* determine the global and max number of subdomains */
2009371c9d4SSatish Balay       struct {
2019371c9d4SSatish Balay         PetscInt max, sum;
2029371c9d4SSatish Balay       } inwork, outwork;
203c10200c1SHong Zhang       PetscMPIInt size;
204c10200c1SHong Zhang 
2056ac3741eSJed Brown       inwork.max = osm->n_local_true;
2066ac3741eSJed Brown       inwork.sum = osm->n_local_true;
207462c564dSBarry Smith       PetscCallMPI(MPIU_Allreduce(&inwork, &outwork, 1, MPIU_2INT, MPIU_MAXSUM_OP, PetscObjectComm((PetscObject)pc)));
2086ac3741eSJed Brown       osm->n_local = outwork.max;
2096ac3741eSJed Brown       osm->n       = outwork.sum;
210c10200c1SHong Zhang 
2119566063dSJacob Faibussowitsch       PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size));
212c10200c1SHong Zhang       if (outwork.max == 1 && outwork.sum == size) {
2137dae84e0SHong Zhang         /* osm->n_local_true = 1 on all processes, set this option may enable use of optimized MatCreateSubMatrices() implementation */
2149566063dSJacob Faibussowitsch         PetscCall(MatSetOption(pc->pmat, MAT_SUBMAT_SINGLEIS, PETSC_TRUE));
215c10200c1SHong Zhang       }
2164b9ad928SBarry Smith     }
21778904715SLisandro Dalcin     if (!osm->is) { /* create the index sets */
2189566063dSJacob Faibussowitsch       PetscCall(PCASMCreateSubdomains(pc->pmat, osm->n_local_true, &osm->is));
2194b9ad928SBarry Smith     }
220f5234e35SJed Brown     if (osm->n_local_true > 1 && !osm->is_local) {
2219566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(osm->n_local_true, &osm->is_local));
222f5234e35SJed Brown       for (i = 0; i < osm->n_local_true; i++) {
223f5234e35SJed Brown         if (osm->overlap > 0) { /* With positive overlap, osm->is[i] will be modified */
2249566063dSJacob Faibussowitsch           PetscCall(ISDuplicate(osm->is[i], &osm->is_local[i]));
2259566063dSJacob Faibussowitsch           PetscCall(ISCopy(osm->is[i], osm->is_local[i]));
226f5234e35SJed Brown         } else {
2279566063dSJacob Faibussowitsch           PetscCall(PetscObjectReference((PetscObject)osm->is[i]));
228f5234e35SJed Brown           osm->is_local[i] = osm->is[i];
229f5234e35SJed Brown         }
230f5234e35SJed Brown       }
231f5234e35SJed Brown     }
2329566063dSJacob Faibussowitsch     PetscCall(PCGetOptionsPrefix(pc, &prefix));
2333d03bb48SJed Brown     if (osm->overlap > 0) {
2344b9ad928SBarry Smith       /* Extend the "overlapping" regions by a number of steps */
2359566063dSJacob Faibussowitsch       PetscCall(MatIncreaseOverlap(pc->pmat, osm->n_local_true, osm->is, osm->overlap));
2363d03bb48SJed Brown     }
2376ed231c7SMatthew Knepley     if (osm->sort_indices) {
23892bb6962SLisandro Dalcin       for (i = 0; i < osm->n_local_true; i++) {
2399566063dSJacob Faibussowitsch         PetscCall(ISSort(osm->is[i]));
24048a46eb9SPierre Jolivet         if (osm->is_local) PetscCall(ISSort(osm->is_local[i]));
2414b9ad928SBarry Smith       }
2426ed231c7SMatthew Knepley     }
24398d3e228SPierre Jolivet     flg = PETSC_FALSE;
24498d3e228SPierre Jolivet     PetscCall(PetscOptionsHasName(NULL, prefix, "-pc_asm_print_subdomains", &flg));
24598d3e228SPierre Jolivet     if (flg) PetscCall(PCASMPrintSubdomains(pc));
246f6991133SBarry Smith     if (!osm->ksp) {
24778904715SLisandro Dalcin       /* Create the local solvers */
2489566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(osm->n_local_true, &osm->ksp));
24948a46eb9SPierre Jolivet       if (domain_dm) PetscCall(PetscInfo(pc, "Setting up ASM subproblems using the embedded DM\n"));
25092bb6962SLisandro Dalcin       for (i = 0; i < osm->n_local_true; i++) {
2519566063dSJacob Faibussowitsch         PetscCall(KSPCreate(PETSC_COMM_SELF, &ksp));
2523821be0aSBarry Smith         PetscCall(KSPSetNestLevel(ksp, pc->kspnestlevel));
2539566063dSJacob Faibussowitsch         PetscCall(KSPSetErrorIfNotConverged(ksp, pc->erroriffailure));
2549566063dSJacob Faibussowitsch         PetscCall(PetscObjectIncrementTabLevel((PetscObject)ksp, (PetscObject)pc, 1));
2559566063dSJacob Faibussowitsch         PetscCall(KSPSetType(ksp, KSPPREONLY));
2569566063dSJacob Faibussowitsch         PetscCall(KSPGetPC(ksp, &subpc));
2579566063dSJacob Faibussowitsch         PetscCall(PCGetOptionsPrefix(pc, &prefix));
2589566063dSJacob Faibussowitsch         PetscCall(KSPSetOptionsPrefix(ksp, prefix));
2599566063dSJacob Faibussowitsch         PetscCall(KSPAppendOptionsPrefix(ksp, "sub_"));
260feb221f8SDmitry Karpeev         if (domain_dm) {
2619566063dSJacob Faibussowitsch           PetscCall(KSPSetDM(ksp, domain_dm[i]));
2629566063dSJacob Faibussowitsch           PetscCall(KSPSetDMActive(ksp, PETSC_FALSE));
2639566063dSJacob Faibussowitsch           PetscCall(DMDestroy(&domain_dm[i]));
264feb221f8SDmitry Karpeev         }
2654b9ad928SBarry Smith         osm->ksp[i] = ksp;
2664b9ad928SBarry Smith       }
2671baa6e33SBarry Smith       if (domain_dm) PetscCall(PetscFree(domain_dm));
268f6991133SBarry Smith     }
2691dd8081eSeaulisa 
2709566063dSJacob Faibussowitsch     PetscCall(ISConcatenate(PETSC_COMM_SELF, osm->n_local_true, osm->is, &osm->lis));
2719566063dSJacob Faibussowitsch     PetscCall(ISSortRemoveDups(osm->lis));
2729566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(osm->lis, &m));
2731dd8081eSeaulisa 
2744b9ad928SBarry Smith     scall = MAT_INITIAL_MATRIX;
2754b9ad928SBarry Smith   } else {
2764b9ad928SBarry Smith     /*
2774b9ad928SBarry Smith        Destroy the blocks from the previous iteration
2784b9ad928SBarry Smith     */
279e09e08ccSBarry Smith     if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
2804849c82aSBarry Smith       PetscCall(MatGetNullSpaces(osm->n_local_true, osm->pmat, &nullsp));
2819566063dSJacob Faibussowitsch       PetscCall(MatDestroyMatrices(osm->n_local_true, &osm->pmat));
2824b9ad928SBarry Smith       scall = MAT_INITIAL_MATRIX;
2834b9ad928SBarry Smith     }
2844b9ad928SBarry Smith   }
2854b9ad928SBarry Smith 
286b58cb649SBarry Smith   /* Destroy previous submatrices of a different type than pc->pmat since MAT_REUSE_MATRIX won't work in that case */
28782b5ce2aSStefano Zampini   if (scall == MAT_REUSE_MATRIX && osm->sub_mat_type) {
2884849c82aSBarry Smith     PetscCall(MatGetNullSpaces(osm->n_local_true, osm->pmat, &nullsp));
28948a46eb9SPierre Jolivet     if (osm->n_local_true > 0) PetscCall(MatDestroySubMatrices(osm->n_local_true, &osm->pmat));
290b58cb649SBarry Smith     scall = MAT_INITIAL_MATRIX;
291b58cb649SBarry Smith   }
292b58cb649SBarry Smith 
29392bb6962SLisandro Dalcin   /*
29492bb6962SLisandro Dalcin      Extract out the submatrices
29592bb6962SLisandro Dalcin   */
2969566063dSJacob Faibussowitsch   PetscCall(MatCreateSubMatrices(pc->pmat, osm->n_local_true, osm->is, osm->is, scall, &osm->pmat));
29792bb6962SLisandro Dalcin   if (scall == MAT_INITIAL_MATRIX) {
2989566063dSJacob Faibussowitsch     PetscCall(PetscObjectGetOptionsPrefix((PetscObject)pc->pmat, &pprefix));
299aa624791SPierre Jolivet     for (i = 0; i < osm->n_local_true; i++) PetscCall(PetscObjectSetOptionsPrefix((PetscObject)osm->pmat[i], pprefix));
3004849c82aSBarry Smith     if (nullsp) PetscCall(MatRestoreNullSpaces(osm->n_local_true, osm->pmat, &nullsp));
30192bb6962SLisandro Dalcin   }
30280ec0b7dSPatrick Sanan 
30380ec0b7dSPatrick Sanan   /* Convert the types of the submatrices (if needbe) */
30480ec0b7dSPatrick Sanan   if (osm->sub_mat_type) {
305f4f49eeaSPierre 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]));
30680ec0b7dSPatrick Sanan   }
30780ec0b7dSPatrick Sanan 
30880ec0b7dSPatrick Sanan   if (!pc->setupcalled) {
30956ea09ceSStefano Zampini     VecType vtype;
31056ea09ceSStefano Zampini 
31180ec0b7dSPatrick Sanan     /* Create the local work vectors (from the local matrices) and scatter contexts */
3129566063dSJacob Faibussowitsch     PetscCall(MatCreateVecs(pc->pmat, &vec, NULL));
3135b3c0d42Seaulisa 
314677c7726SPierre Jolivet     PetscCheck(!osm->is_local || osm->n_local_true == 1 || (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() with more than a single subdomain");
315677c7726SPierre Jolivet     if (osm->is_local && osm->type != PC_ASM_BASIC && osm->loctype == PC_COMPOSITE_ADDITIVE) PetscCall(PetscMalloc1(osm->n_local_true, &osm->lprolongation));
3169566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(osm->n_local_true, &osm->lrestriction));
3179566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(osm->n_local_true, &osm->x));
3189566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(osm->n_local_true, &osm->y));
319347574c9Seaulisa 
3209566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(osm->lis, &m));
3219566063dSJacob Faibussowitsch     PetscCall(ISCreateStride(PETSC_COMM_SELF, m, 0, 1, &isl));
3229566063dSJacob Faibussowitsch     PetscCall(MatGetVecType(osm->pmat[0], &vtype));
3239566063dSJacob Faibussowitsch     PetscCall(VecCreate(PETSC_COMM_SELF, &osm->lx));
3249566063dSJacob Faibussowitsch     PetscCall(VecSetSizes(osm->lx, m, m));
3259566063dSJacob Faibussowitsch     PetscCall(VecSetType(osm->lx, vtype));
3269566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(osm->lx, &osm->ly));
3279566063dSJacob Faibussowitsch     PetscCall(VecScatterCreate(vec, osm->lis, osm->lx, isl, &osm->restriction));
3289566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&isl));
329347574c9Seaulisa 
33080ec0b7dSPatrick Sanan     for (i = 0; i < osm->n_local_true; ++i) {
3315b3c0d42Seaulisa       ISLocalToGlobalMapping ltog;
3325b3c0d42Seaulisa       IS                     isll;
3335b3c0d42Seaulisa       const PetscInt        *idx_is;
3345b3c0d42Seaulisa       PetscInt              *idx_lis, nout;
3355b3c0d42Seaulisa 
3369566063dSJacob Faibussowitsch       PetscCall(ISGetLocalSize(osm->is[i], &m));
3379566063dSJacob Faibussowitsch       PetscCall(MatCreateVecs(osm->pmat[i], &osm->x[i], NULL));
3389566063dSJacob Faibussowitsch       PetscCall(VecDuplicate(osm->x[i], &osm->y[i]));
33955817e79SBarry Smith 
340b0de9f23SBarry Smith       /* generate a scatter from ly to y[i] picking all the overlapping is[i] entries */
3419566063dSJacob Faibussowitsch       PetscCall(ISLocalToGlobalMappingCreateIS(osm->lis, &ltog));
3429566063dSJacob Faibussowitsch       PetscCall(ISGetLocalSize(osm->is[i], &m));
3439566063dSJacob Faibussowitsch       PetscCall(ISGetIndices(osm->is[i], &idx_is));
3449566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(m, &idx_lis));
3459566063dSJacob Faibussowitsch       PetscCall(ISGlobalToLocalMappingApply(ltog, IS_GTOLM_DROP, m, idx_is, &nout, idx_lis));
3467827d75bSBarry Smith       PetscCheck(nout == m, PETSC_COMM_SELF, PETSC_ERR_PLIB, "is not a subset of lis");
3479566063dSJacob Faibussowitsch       PetscCall(ISRestoreIndices(osm->is[i], &idx_is));
3489566063dSJacob Faibussowitsch       PetscCall(ISCreateGeneral(PETSC_COMM_SELF, m, idx_lis, PETSC_OWN_POINTER, &isll));
3499566063dSJacob Faibussowitsch       PetscCall(ISLocalToGlobalMappingDestroy(&ltog));
3509566063dSJacob Faibussowitsch       PetscCall(ISCreateStride(PETSC_COMM_SELF, m, 0, 1, &isl));
3519566063dSJacob Faibussowitsch       PetscCall(VecScatterCreate(osm->ly, isll, osm->y[i], isl, &osm->lrestriction[i]));
3529566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&isll));
3539566063dSJacob Faibussowitsch       PetscCall(ISDestroy(&isl));
35415229ffcSPierre Jolivet       if (osm->lprolongation) { /* generate a scatter from y[i] to ly picking only the non-overlapping is_local[i] entries */
355d8b3b5e3Seaulisa         ISLocalToGlobalMapping ltog;
356d8b3b5e3Seaulisa         IS                     isll, isll_local;
357d8b3b5e3Seaulisa         const PetscInt        *idx_local;
358d8b3b5e3Seaulisa         PetscInt              *idx1, *idx2, nout;
359d8b3b5e3Seaulisa 
3609566063dSJacob Faibussowitsch         PetscCall(ISGetLocalSize(osm->is_local[i], &m_local));
3619566063dSJacob Faibussowitsch         PetscCall(ISGetIndices(osm->is_local[i], &idx_local));
362d8b3b5e3Seaulisa 
3639566063dSJacob Faibussowitsch         PetscCall(ISLocalToGlobalMappingCreateIS(osm->is[i], &ltog));
3649566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(m_local, &idx1));
3659566063dSJacob Faibussowitsch         PetscCall(ISGlobalToLocalMappingApply(ltog, IS_GTOLM_DROP, m_local, idx_local, &nout, idx1));
3669566063dSJacob Faibussowitsch         PetscCall(ISLocalToGlobalMappingDestroy(&ltog));
3677827d75bSBarry Smith         PetscCheck(nout == m_local, PETSC_COMM_SELF, PETSC_ERR_PLIB, "is_local not a subset of is");
3689566063dSJacob Faibussowitsch         PetscCall(ISCreateGeneral(PETSC_COMM_SELF, m_local, idx1, PETSC_OWN_POINTER, &isll));
369d8b3b5e3Seaulisa 
3709566063dSJacob Faibussowitsch         PetscCall(ISLocalToGlobalMappingCreateIS(osm->lis, &ltog));
3719566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(m_local, &idx2));
3729566063dSJacob Faibussowitsch         PetscCall(ISGlobalToLocalMappingApply(ltog, IS_GTOLM_DROP, m_local, idx_local, &nout, idx2));
3739566063dSJacob Faibussowitsch         PetscCall(ISLocalToGlobalMappingDestroy(&ltog));
3747827d75bSBarry Smith         PetscCheck(nout == m_local, PETSC_COMM_SELF, PETSC_ERR_PLIB, "is_local not a subset of lis");
3759566063dSJacob Faibussowitsch         PetscCall(ISCreateGeneral(PETSC_COMM_SELF, m_local, idx2, PETSC_OWN_POINTER, &isll_local));
376d8b3b5e3Seaulisa 
3779566063dSJacob Faibussowitsch         PetscCall(ISRestoreIndices(osm->is_local[i], &idx_local));
3789566063dSJacob Faibussowitsch         PetscCall(VecScatterCreate(osm->y[i], isll, osm->ly, isll_local, &osm->lprolongation[i]));
379d8b3b5e3Seaulisa 
3809566063dSJacob Faibussowitsch         PetscCall(ISDestroy(&isll));
3819566063dSJacob Faibussowitsch         PetscCall(ISDestroy(&isll_local));
382d8b3b5e3Seaulisa       }
38380ec0b7dSPatrick Sanan     }
3849566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&vec));
38580ec0b7dSPatrick Sanan   }
38680ec0b7dSPatrick Sanan 
387fb745f2cSMatthew G. Knepley   if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE) {
388235cc4ceSMatthew G. Knepley     IS      *cis;
389235cc4ceSMatthew G. Knepley     PetscInt c;
390235cc4ceSMatthew G. Knepley 
3919566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(osm->n_local_true, &cis));
392235cc4ceSMatthew G. Knepley     for (c = 0; c < osm->n_local_true; ++c) cis[c] = osm->lis;
3939566063dSJacob Faibussowitsch     PetscCall(MatCreateSubMatrices(pc->pmat, osm->n_local_true, osm->is, cis, scall, &osm->lmats));
3949566063dSJacob Faibussowitsch     PetscCall(PetscFree(cis));
395fb745f2cSMatthew G. Knepley   }
3964b9ad928SBarry Smith 
3974b9ad928SBarry Smith   /* Return control to the user so that the submatrices can be modified (e.g., to apply
3984b9ad928SBarry Smith      different boundary conditions for the submatrices than for the global problem) */
3999566063dSJacob Faibussowitsch   PetscCall(PCModifySubMatrices(pc, osm->n_local_true, osm->is, osm->is, osm->pmat, pc->modifysubmatricesP));
4004b9ad928SBarry Smith 
40192bb6962SLisandro Dalcin   /*
40292bb6962SLisandro Dalcin      Loop over subdomains putting them into local ksp
40392bb6962SLisandro Dalcin   */
4049566063dSJacob Faibussowitsch   PetscCall(KSPGetOptionsPrefix(osm->ksp[0], &prefix));
40592bb6962SLisandro Dalcin   for (i = 0; i < osm->n_local_true; i++) {
4069566063dSJacob Faibussowitsch     PetscCall(KSPSetOperators(osm->ksp[i], osm->pmat[i], osm->pmat[i]));
4079566063dSJacob Faibussowitsch     PetscCall(MatSetOptionsPrefix(osm->pmat[i], prefix));
40848a46eb9SPierre Jolivet     if (!pc->setupcalled) PetscCall(KSPSetFromOptions(osm->ksp[i]));
40992bb6962SLisandro Dalcin   }
4103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4114b9ad928SBarry Smith }
4124b9ad928SBarry Smith 
413d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUpOnBlocks_ASM(PC pc)
414d71ae5a4SJacob Faibussowitsch {
4154b9ad928SBarry Smith   PC_ASM            *osm = (PC_ASM *)pc->data;
41613f74950SBarry Smith   PetscInt           i;
417539c167fSBarry Smith   KSPConvergedReason reason;
4184b9ad928SBarry Smith 
4194b9ad928SBarry Smith   PetscFunctionBegin;
4204b9ad928SBarry Smith   for (i = 0; i < osm->n_local_true; i++) {
4219566063dSJacob Faibussowitsch     PetscCall(KSPSetUp(osm->ksp[i]));
4229566063dSJacob Faibussowitsch     PetscCall(KSPGetConvergedReason(osm->ksp[i], &reason));
423ad540459SPierre Jolivet     if (reason == KSP_DIVERGED_PC_FAILED) pc->failedreason = PC_SUBPC_ERROR;
4244b9ad928SBarry Smith   }
4253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4264b9ad928SBarry Smith }
4274b9ad928SBarry Smith 
428d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_ASM(PC pc, Vec x, Vec y)
429d71ae5a4SJacob Faibussowitsch {
4304b9ad928SBarry Smith   PC_ASM     *osm = (PC_ASM *)pc->data;
4311dd8081eSeaulisa   PetscInt    i, n_local_true = osm->n_local_true;
4324b9ad928SBarry Smith   ScatterMode forward = SCATTER_FORWARD, reverse = SCATTER_REVERSE;
4334b9ad928SBarry Smith 
4344b9ad928SBarry Smith   PetscFunctionBegin;
4354b9ad928SBarry Smith   /*
43648e38a8aSPierre Jolivet      support for limiting the restriction or interpolation to only local
4374b9ad928SBarry Smith      subdomain values (leaving the other values 0).
4384b9ad928SBarry Smith   */
4394b9ad928SBarry Smith   if (!(osm->type & PC_ASM_RESTRICT)) {
4404b9ad928SBarry Smith     forward = SCATTER_FORWARD_LOCAL;
4414b9ad928SBarry Smith     /* have to zero the work RHS since scatter may leave some slots empty */
4429566063dSJacob Faibussowitsch     PetscCall(VecSet(osm->lx, 0.0));
4434b9ad928SBarry Smith   }
444ad540459SPierre Jolivet   if (!(osm->type & PC_ASM_INTERPOLATE)) reverse = SCATTER_REVERSE_LOCAL;
4454b9ad928SBarry Smith 
4460fdf79fbSJacob Faibussowitsch   PetscCheck(osm->loctype == PC_COMPOSITE_MULTIPLICATIVE || osm->loctype == PC_COMPOSITE_ADDITIVE, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Invalid local composition type: %s", PCCompositeTypes[osm->loctype]);
447b0de9f23SBarry Smith   /* zero the global and the local solutions */
4489566063dSJacob Faibussowitsch   PetscCall(VecSet(y, 0.0));
4499566063dSJacob Faibussowitsch   PetscCall(VecSet(osm->ly, 0.0));
450347574c9Seaulisa 
45148e38a8aSPierre Jolivet   /* copy the global RHS to local RHS including the ghost nodes */
4529566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(osm->restriction, x, osm->lx, INSERT_VALUES, forward));
4539566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(osm->restriction, x, osm->lx, INSERT_VALUES, forward));
454347574c9Seaulisa 
45548e38a8aSPierre Jolivet   /* restrict local RHS to the overlapping 0-block RHS */
4569566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward));
4579566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward));
458d8b3b5e3Seaulisa 
45912cd4985SMatthew G. Knepley   /* do the local solves */
46012cd4985SMatthew G. Knepley   for (i = 0; i < n_local_true; ++i) {
461b0de9f23SBarry Smith     /* solve the overlapping i-block */
4629566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(PC_ApplyOnBlocks, osm->ksp[i], osm->x[i], osm->y[i], 0));
4639566063dSJacob Faibussowitsch     PetscCall(KSPSolve(osm->ksp[i], osm->x[i], osm->y[i]));
4649566063dSJacob Faibussowitsch     PetscCall(KSPCheckSolve(osm->ksp[i], pc, osm->y[i]));
4659566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(PC_ApplyOnBlocks, osm->ksp[i], osm->x[i], osm->y[i], 0));
466d8b3b5e3Seaulisa 
467677c7726SPierre Jolivet     if (osm->lprolongation && osm->type != PC_ASM_INTERPOLATE) { /* interpolate the non-overlapping i-block solution to the local solution (only for restrictive additive) */
4689566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward));
4699566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward));
47048e38a8aSPierre Jolivet     } else { /* interpolate the overlapping i-block solution to the local solution */
4719566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse));
4729566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse));
473d8b3b5e3Seaulisa     }
474347574c9Seaulisa 
475347574c9Seaulisa     if (i < n_local_true - 1) {
47648e38a8aSPierre Jolivet       /* restrict local RHS to the overlapping (i+1)-block RHS */
4779566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(osm->lrestriction[i + 1], osm->lx, osm->x[i + 1], INSERT_VALUES, forward));
4789566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(osm->lrestriction[i + 1], osm->lx, osm->x[i + 1], INSERT_VALUES, forward));
479347574c9Seaulisa 
480347574c9Seaulisa       if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE) {
4818966356dSPierre Jolivet         /* update the overlapping (i+1)-block RHS using the current local solution */
4829566063dSJacob Faibussowitsch         PetscCall(MatMult(osm->lmats[i + 1], osm->ly, osm->y[i + 1]));
4839566063dSJacob Faibussowitsch         PetscCall(VecAXPBY(osm->x[i + 1], -1., 1., osm->y[i + 1]));
4847c3d802fSMatthew G. Knepley       }
48512cd4985SMatthew G. Knepley     }
48612cd4985SMatthew G. Knepley   }
48748e38a8aSPierre Jolivet   /* add the local solution to the global solution including the ghost nodes */
4889566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(osm->restriction, osm->ly, y, ADD_VALUES, reverse));
4899566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(osm->restriction, osm->ly, y, ADD_VALUES, reverse));
4903ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
4914b9ad928SBarry Smith }
4924b9ad928SBarry Smith 
493*4dbf25a8SPierre Jolivet static PetscErrorCode PCMatApply_ASM_Private(PC pc, Mat X, Mat Y, PetscBool transpose)
494d71ae5a4SJacob Faibussowitsch {
49548e38a8aSPierre Jolivet   PC_ASM     *osm = (PC_ASM *)pc->data;
49648e38a8aSPierre Jolivet   Mat         Z, W;
49748e38a8aSPierre Jolivet   Vec         x;
49848e38a8aSPierre Jolivet   PetscInt    i, m, N;
49948e38a8aSPierre Jolivet   ScatterMode forward = SCATTER_FORWARD, reverse = SCATTER_REVERSE;
50048e38a8aSPierre Jolivet 
50148e38a8aSPierre Jolivet   PetscFunctionBegin;
5027827d75bSBarry Smith   PetscCheck(osm->n_local_true <= 1, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Not yet implemented");
50348e38a8aSPierre Jolivet   /*
50448e38a8aSPierre Jolivet      support for limiting the restriction or interpolation to only local
50548e38a8aSPierre Jolivet      subdomain values (leaving the other values 0).
50648e38a8aSPierre Jolivet   */
507*4dbf25a8SPierre Jolivet   if ((!transpose && !(osm->type & PC_ASM_RESTRICT)) || (transpose && !(osm->type & PC_ASM_INTERPOLATE))) {
50848e38a8aSPierre Jolivet     forward = SCATTER_FORWARD_LOCAL;
50948e38a8aSPierre Jolivet     /* have to zero the work RHS since scatter may leave some slots empty */
5109566063dSJacob Faibussowitsch     PetscCall(VecSet(osm->lx, 0.0));
51148e38a8aSPierre Jolivet   }
512*4dbf25a8SPierre Jolivet   if ((!transpose && !(osm->type & PC_ASM_INTERPOLATE)) || (transpose && !(osm->type & PC_ASM_RESTRICT))) reverse = SCATTER_REVERSE_LOCAL;
5139566063dSJacob Faibussowitsch   PetscCall(VecGetLocalSize(osm->x[0], &m));
5149566063dSJacob Faibussowitsch   PetscCall(MatGetSize(X, NULL, &N));
5159566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, m, N, NULL, &Z));
5160fdf79fbSJacob Faibussowitsch 
5170fdf79fbSJacob Faibussowitsch   PetscCheck(osm->loctype == PC_COMPOSITE_MULTIPLICATIVE || osm->loctype == PC_COMPOSITE_ADDITIVE, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Invalid local composition type: %s", PCCompositeTypes[osm->loctype]);
51848e38a8aSPierre Jolivet   /* zero the global and the local solutions */
5199566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(Y));
5209566063dSJacob Faibussowitsch   PetscCall(VecSet(osm->ly, 0.0));
52148e38a8aSPierre Jolivet 
52248e38a8aSPierre Jolivet   for (i = 0; i < N; ++i) {
5239566063dSJacob Faibussowitsch     PetscCall(MatDenseGetColumnVecRead(X, i, &x));
52448e38a8aSPierre Jolivet     /* copy the global RHS to local RHS including the ghost nodes */
5259566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(osm->restriction, x, osm->lx, INSERT_VALUES, forward));
5269566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(osm->restriction, x, osm->lx, INSERT_VALUES, forward));
5279566063dSJacob Faibussowitsch     PetscCall(MatDenseRestoreColumnVecRead(X, i, &x));
52848e38a8aSPierre Jolivet 
5299566063dSJacob Faibussowitsch     PetscCall(MatDenseGetColumnVecWrite(Z, i, &x));
53048e38a8aSPierre Jolivet     /* restrict local RHS to the overlapping 0-block RHS */
5319566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(osm->lrestriction[0], osm->lx, x, INSERT_VALUES, forward));
5329566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(osm->lrestriction[0], osm->lx, x, INSERT_VALUES, forward));
5339566063dSJacob Faibussowitsch     PetscCall(MatDenseRestoreColumnVecWrite(Z, i, &x));
53448e38a8aSPierre Jolivet   }
5359566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, m, N, NULL, &W));
53648e38a8aSPierre Jolivet   /* solve the overlapping 0-block */
537*4dbf25a8SPierre Jolivet   if (!transpose) {
5389566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(PC_ApplyOnBlocks, osm->ksp[0], Z, W, 0));
5399566063dSJacob Faibussowitsch     PetscCall(KSPMatSolve(osm->ksp[0], Z, W));
5409566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(PC_ApplyOnBlocks, osm->ksp[0], Z, W, 0));
541*4dbf25a8SPierre Jolivet   } else {
542*4dbf25a8SPierre Jolivet     PetscCall(PetscLogEventBegin(PC_ApplyTransposeOnBlocks, osm->ksp[0], Z, W, 0));
543*4dbf25a8SPierre Jolivet     PetscCall(KSPMatSolveTranspose(osm->ksp[0], Z, W));
544*4dbf25a8SPierre Jolivet     PetscCall(PetscLogEventEnd(PC_ApplyTransposeOnBlocks, osm->ksp[0], Z, W, 0));
545*4dbf25a8SPierre Jolivet   }
546*4dbf25a8SPierre Jolivet   PetscCall(KSPCheckSolve(osm->ksp[0], pc, NULL));
5479566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Z));
54848e38a8aSPierre Jolivet 
54948e38a8aSPierre Jolivet   for (i = 0; i < N; ++i) {
5509566063dSJacob Faibussowitsch     PetscCall(VecSet(osm->ly, 0.0));
5519566063dSJacob Faibussowitsch     PetscCall(MatDenseGetColumnVecRead(W, i, &x));
552*4dbf25a8SPierre Jolivet     if (osm->lprolongation && ((!transpose && osm->type != PC_ASM_INTERPOLATE) || (transpose && osm->type != PC_ASM_RESTRICT))) { /* interpolate the non-overlapping 0-block solution to the local solution (only for restrictive additive) */
5539566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(osm->lprolongation[0], x, osm->ly, ADD_VALUES, forward));
5549566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(osm->lprolongation[0], x, osm->ly, ADD_VALUES, forward));
55548e38a8aSPierre Jolivet     } else { /* interpolate the overlapping 0-block solution to the local solution */
5569566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(osm->lrestriction[0], x, osm->ly, ADD_VALUES, reverse));
5579566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(osm->lrestriction[0], x, osm->ly, ADD_VALUES, reverse));
55848e38a8aSPierre Jolivet     }
5599566063dSJacob Faibussowitsch     PetscCall(MatDenseRestoreColumnVecRead(W, i, &x));
56048e38a8aSPierre Jolivet 
5619566063dSJacob Faibussowitsch     PetscCall(MatDenseGetColumnVecWrite(Y, i, &x));
56248e38a8aSPierre Jolivet     /* add the local solution to the global solution including the ghost nodes */
5639566063dSJacob Faibussowitsch     PetscCall(VecScatterBegin(osm->restriction, osm->ly, x, ADD_VALUES, reverse));
5649566063dSJacob Faibussowitsch     PetscCall(VecScatterEnd(osm->restriction, osm->ly, x, ADD_VALUES, reverse));
5659566063dSJacob Faibussowitsch     PetscCall(MatDenseRestoreColumnVecWrite(Y, i, &x));
56648e38a8aSPierre Jolivet   }
5679566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&W));
5683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
56948e38a8aSPierre Jolivet }
57048e38a8aSPierre Jolivet 
571*4dbf25a8SPierre Jolivet static PetscErrorCode PCMatApply_ASM(PC pc, Mat X, Mat Y)
572*4dbf25a8SPierre Jolivet {
573*4dbf25a8SPierre Jolivet   PetscFunctionBegin;
574*4dbf25a8SPierre Jolivet   PetscCall(PCMatApply_ASM_Private(pc, X, Y, PETSC_FALSE));
575*4dbf25a8SPierre Jolivet   PetscFunctionReturn(PETSC_SUCCESS);
576*4dbf25a8SPierre Jolivet }
577*4dbf25a8SPierre Jolivet 
578*4dbf25a8SPierre Jolivet static PetscErrorCode PCMatApplyTranspose_ASM(PC pc, Mat X, Mat Y)
579*4dbf25a8SPierre Jolivet {
580*4dbf25a8SPierre Jolivet   PetscFunctionBegin;
581*4dbf25a8SPierre Jolivet   PetscCall(PCMatApply_ASM_Private(pc, X, Y, PETSC_TRUE));
582*4dbf25a8SPierre Jolivet   PetscFunctionReturn(PETSC_SUCCESS);
583*4dbf25a8SPierre Jolivet }
584*4dbf25a8SPierre Jolivet 
585d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplyTranspose_ASM(PC pc, Vec x, Vec y)
586d71ae5a4SJacob Faibussowitsch {
5874b9ad928SBarry Smith   PC_ASM     *osm = (PC_ASM *)pc->data;
5881dd8081eSeaulisa   PetscInt    i, n_local_true = osm->n_local_true;
5894b9ad928SBarry Smith   ScatterMode forward = SCATTER_FORWARD, reverse = SCATTER_REVERSE;
5904b9ad928SBarry Smith 
5914b9ad928SBarry Smith   PetscFunctionBegin;
5924b9ad928SBarry Smith   /*
5934b9ad928SBarry Smith      Support for limiting the restriction or interpolation to only local
5944b9ad928SBarry Smith      subdomain values (leaving the other values 0).
5954b9ad928SBarry Smith 
5964b9ad928SBarry Smith      Note: these are reversed from the PCApply_ASM() because we are applying the
5974b9ad928SBarry Smith      transpose of the three terms
5984b9ad928SBarry Smith   */
599d8b3b5e3Seaulisa 
6004b9ad928SBarry Smith   if (!(osm->type & PC_ASM_INTERPOLATE)) {
6014b9ad928SBarry Smith     forward = SCATTER_FORWARD_LOCAL;
6024b9ad928SBarry Smith     /* have to zero the work RHS since scatter may leave some slots empty */
6039566063dSJacob Faibussowitsch     PetscCall(VecSet(osm->lx, 0.0));
6044b9ad928SBarry Smith   }
6052fa5cd67SKarl Rupp   if (!(osm->type & PC_ASM_RESTRICT)) reverse = SCATTER_REVERSE_LOCAL;
6064b9ad928SBarry Smith 
607b0de9f23SBarry Smith   /* zero the global and the local solutions */
6089566063dSJacob Faibussowitsch   PetscCall(VecSet(y, 0.0));
6099566063dSJacob Faibussowitsch   PetscCall(VecSet(osm->ly, 0.0));
610d8b3b5e3Seaulisa 
611b0de9f23SBarry Smith   /* Copy the global RHS to local RHS including the ghost nodes */
6129566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(osm->restriction, x, osm->lx, INSERT_VALUES, forward));
6139566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(osm->restriction, x, osm->lx, INSERT_VALUES, forward));
614d8b3b5e3Seaulisa 
615b0de9f23SBarry Smith   /* Restrict local RHS to the overlapping 0-block RHS */
6169566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward));
6179566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward));
618d8b3b5e3Seaulisa 
6194b9ad928SBarry Smith   /* do the local solves */
620d8b3b5e3Seaulisa   for (i = 0; i < n_local_true; ++i) {
621b0de9f23SBarry Smith     /* solve the overlapping i-block */
6229566063dSJacob Faibussowitsch     PetscCall(PetscLogEventBegin(PC_ApplyOnBlocks, osm->ksp[i], osm->x[i], osm->y[i], 0));
6239566063dSJacob Faibussowitsch     PetscCall(KSPSolveTranspose(osm->ksp[i], osm->x[i], osm->y[i]));
6249566063dSJacob Faibussowitsch     PetscCall(KSPCheckSolve(osm->ksp[i], pc, osm->y[i]));
6259566063dSJacob Faibussowitsch     PetscCall(PetscLogEventEnd(PC_ApplyOnBlocks, osm->ksp[i], osm->x[i], osm->y[i], 0));
626d8b3b5e3Seaulisa 
627a681a0f1SPierre Jolivet     if (osm->lprolongation && osm->type != PC_ASM_RESTRICT) { /* interpolate the non-overlapping i-block solution to the local solution */
6289566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward));
6299566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward));
630910cf402Sprj-     } else { /* interpolate the overlapping i-block solution to the local solution */
6319566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse));
6329566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse));
6334b9ad928SBarry Smith     }
634d8b3b5e3Seaulisa 
635d8b3b5e3Seaulisa     if (i < n_local_true - 1) {
636b0de9f23SBarry Smith       /* Restrict local RHS to the overlapping (i+1)-block RHS */
6379566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(osm->lrestriction[i + 1], osm->lx, osm->x[i + 1], INSERT_VALUES, forward));
6389566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(osm->lrestriction[i + 1], osm->lx, osm->x[i + 1], INSERT_VALUES, forward));
6394b9ad928SBarry Smith     }
6404b9ad928SBarry Smith   }
641b0de9f23SBarry Smith   /* Add the local solution to the global solution including the ghost nodes */
6429566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(osm->restriction, osm->ly, y, ADD_VALUES, reverse));
6439566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(osm->restriction, osm->ly, y, ADD_VALUES, reverse));
6443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6454b9ad928SBarry Smith }
6464b9ad928SBarry Smith 
647d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCReset_ASM(PC pc)
648d71ae5a4SJacob Faibussowitsch {
6494b9ad928SBarry Smith   PC_ASM  *osm = (PC_ASM *)pc->data;
65013f74950SBarry Smith   PetscInt i;
6514b9ad928SBarry Smith 
6524b9ad928SBarry Smith   PetscFunctionBegin;
65392bb6962SLisandro Dalcin   if (osm->ksp) {
65448a46eb9SPierre Jolivet     for (i = 0; i < osm->n_local_true; i++) PetscCall(KSPReset(osm->ksp[i]));
65592bb6962SLisandro Dalcin   }
656e09e08ccSBarry Smith   if (osm->pmat) {
65748a46eb9SPierre Jolivet     if (osm->n_local_true > 0) PetscCall(MatDestroySubMatrices(osm->n_local_true, &osm->pmat));
65892bb6962SLisandro Dalcin   }
6591dd8081eSeaulisa   if (osm->lrestriction) {
6609566063dSJacob Faibussowitsch     PetscCall(VecScatterDestroy(&osm->restriction));
6611dd8081eSeaulisa     for (i = 0; i < osm->n_local_true; i++) {
6629566063dSJacob Faibussowitsch       PetscCall(VecScatterDestroy(&osm->lrestriction[i]));
6639566063dSJacob Faibussowitsch       if (osm->lprolongation) PetscCall(VecScatterDestroy(&osm->lprolongation[i]));
6649566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&osm->x[i]));
6659566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&osm->y[i]));
6664b9ad928SBarry Smith     }
6679566063dSJacob Faibussowitsch     PetscCall(PetscFree(osm->lrestriction));
6689566063dSJacob Faibussowitsch     if (osm->lprolongation) PetscCall(PetscFree(osm->lprolongation));
6699566063dSJacob Faibussowitsch     PetscCall(PetscFree(osm->x));
6709566063dSJacob Faibussowitsch     PetscCall(PetscFree(osm->y));
67192bb6962SLisandro Dalcin   }
672ce78bad3SBarry Smith   PetscCall(PCASMDestroySubdomains(osm->n_local_true, &osm->is, &osm->is_local));
6739566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&osm->lis));
6749566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&osm->lx));
6759566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&osm->ly));
67648a46eb9SPierre Jolivet   if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE) PetscCall(MatDestroyMatrices(osm->n_local_true, &osm->lmats));
6772fa5cd67SKarl Rupp 
6789566063dSJacob Faibussowitsch   PetscCall(PetscFree(osm->sub_mat_type));
67980ec0b7dSPatrick Sanan 
6800a545947SLisandro Dalcin   osm->is       = NULL;
6810a545947SLisandro Dalcin   osm->is_local = NULL;
6823ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
683e91c6855SBarry Smith }
684e91c6855SBarry Smith 
685d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_ASM(PC pc)
686d71ae5a4SJacob Faibussowitsch {
687e91c6855SBarry Smith   PC_ASM  *osm = (PC_ASM *)pc->data;
688e91c6855SBarry Smith   PetscInt i;
689e91c6855SBarry Smith 
690e91c6855SBarry Smith   PetscFunctionBegin;
6919566063dSJacob Faibussowitsch   PetscCall(PCReset_ASM(pc));
692e91c6855SBarry Smith   if (osm->ksp) {
69348a46eb9SPierre Jolivet     for (i = 0; i < osm->n_local_true; i++) PetscCall(KSPDestroy(&osm->ksp[i]));
6949566063dSJacob Faibussowitsch     PetscCall(PetscFree(osm->ksp));
695e91c6855SBarry Smith   }
6969566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc->data));
69796322394SPierre Jolivet 
6989566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetLocalSubdomains_C", NULL));
6999566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetTotalSubdomains_C", NULL));
7009566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetOverlap_C", NULL));
7019566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetType_C", NULL));
7029566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMGetType_C", NULL));
7039566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetLocalType_C", NULL));
7049566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMGetLocalType_C", NULL));
7059566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetSortIndices_C", NULL));
7069566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMGetSubKSP_C", NULL));
7079566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMGetSubMatType_C", NULL));
7089566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetSubMatType_C", NULL));
7093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7104b9ad928SBarry Smith }
7114b9ad928SBarry Smith 
712ce78bad3SBarry Smith static PetscErrorCode PCSetFromOptions_ASM(PC pc, PetscOptionItems PetscOptionsObject)
713d71ae5a4SJacob Faibussowitsch {
7144b9ad928SBarry Smith   PC_ASM         *osm = (PC_ASM *)pc->data;
7159dcbbd2bSBarry Smith   PetscInt        blocks, ovl;
71657501b6eSBarry Smith   PetscBool       flg;
71792bb6962SLisandro Dalcin   PCASMType       asmtype;
71812cd4985SMatthew G. Knepley   PCCompositeType loctype;
71980ec0b7dSPatrick Sanan   char            sub_mat_type[256];
7204b9ad928SBarry Smith 
7214b9ad928SBarry Smith   PetscFunctionBegin;
722d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Additive Schwarz options");
7239566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_asm_dm_subdomains", "Use DMCreateDomainDecomposition() to define subdomains", "PCASMSetDMSubdomains", osm->dm_subdomains, &osm->dm_subdomains, &flg));
7249566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_asm_blocks", "Number of subdomains", "PCASMSetTotalSubdomains", osm->n, &blocks, &flg));
72565db9045SDmitry Karpeev   if (flg) {
7269566063dSJacob Faibussowitsch     PetscCall(PCASMSetTotalSubdomains(pc, blocks, NULL, NULL));
727d709ea83SDmitry Karpeev     osm->dm_subdomains = PETSC_FALSE;
72865db9045SDmitry Karpeev   }
7299566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_asm_local_blocks", "Number of local subdomains", "PCASMSetLocalSubdomains", osm->n_local_true, &blocks, &flg));
730342c94f9SMatthew G. Knepley   if (flg) {
7319566063dSJacob Faibussowitsch     PetscCall(PCASMSetLocalSubdomains(pc, blocks, NULL, NULL));
732342c94f9SMatthew G. Knepley     osm->dm_subdomains = PETSC_FALSE;
733342c94f9SMatthew G. Knepley   }
7349566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_asm_overlap", "Number of grid points overlap", "PCASMSetOverlap", osm->overlap, &ovl, &flg));
73565db9045SDmitry Karpeev   if (flg) {
7369566063dSJacob Faibussowitsch     PetscCall(PCASMSetOverlap(pc, ovl));
737d709ea83SDmitry Karpeev     osm->dm_subdomains = PETSC_FALSE;
73865db9045SDmitry Karpeev   }
73990d69ab7SBarry Smith   flg = PETSC_FALSE;
7409566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-pc_asm_type", "Type of restriction/extension", "PCASMSetType", PCASMTypes, (PetscEnum)osm->type, (PetscEnum *)&asmtype, &flg));
7419566063dSJacob Faibussowitsch   if (flg) PetscCall(PCASMSetType(pc, asmtype));
74212cd4985SMatthew G. Knepley   flg = PETSC_FALSE;
7439566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-pc_asm_local_type", "Type of local solver composition", "PCASMSetLocalType", PCCompositeTypes, (PetscEnum)osm->loctype, (PetscEnum *)&loctype, &flg));
7449566063dSJacob Faibussowitsch   if (flg) PetscCall(PCASMSetLocalType(pc, loctype));
7459566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFList("-pc_asm_sub_mat_type", "Subsolve Matrix Type", "PCASMSetSubMatType", MatList, NULL, sub_mat_type, 256, &flg));
7461baa6e33SBarry Smith   if (flg) PetscCall(PCASMSetSubMatType(pc, sub_mat_type));
747d0609cedSBarry Smith   PetscOptionsHeadEnd();
7483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7494b9ad928SBarry Smith }
7504b9ad928SBarry Smith 
751d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCASMSetLocalSubdomains_ASM(PC pc, PetscInt n, IS is[], IS is_local[])
752d71ae5a4SJacob Faibussowitsch {
7534b9ad928SBarry Smith   PC_ASM  *osm = (PC_ASM *)pc->data;
75492bb6962SLisandro Dalcin   PetscInt i;
7554b9ad928SBarry Smith 
7564b9ad928SBarry Smith   PetscFunctionBegin;
75763a3b9bcSJacob Faibussowitsch   PetscCheck(n >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Each process must have 1 or more blocks, n = %" PetscInt_FMT, n);
7587827d75bSBarry Smith   PetscCheck(!pc->setupcalled || (n == osm->n_local_true && !is), PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "PCASMSetLocalSubdomains() should be called before calling PCSetUp().");
759e7e72b3dSBarry Smith 
7604b9ad928SBarry Smith   if (!pc->setupcalled) {
76192bb6962SLisandro Dalcin     if (is) {
7629566063dSJacob Faibussowitsch       for (i = 0; i < n; i++) PetscCall(PetscObjectReference((PetscObject)is[i]));
76392bb6962SLisandro Dalcin     }
764832fc9a5SMatthew Knepley     if (is_local) {
7659566063dSJacob Faibussowitsch       for (i = 0; i < n; i++) PetscCall(PetscObjectReference((PetscObject)is_local[i]));
766832fc9a5SMatthew Knepley     }
767ce78bad3SBarry Smith     PetscCall(PCASMDestroySubdomains(osm->n_local_true, &osm->is, &osm->is_local));
7682fa5cd67SKarl Rupp 
76907517c86SMark Adams     if (osm->ksp && osm->n_local_true != n) {
77007517c86SMark Adams       for (i = 0; i < osm->n_local_true; i++) PetscCall(KSPDestroy(&osm->ksp[i]));
77107517c86SMark Adams       PetscCall(PetscFree(osm->ksp));
77207517c86SMark Adams     }
77307517c86SMark Adams 
7744b9ad928SBarry Smith     osm->n_local_true = n;
7750a545947SLisandro Dalcin     osm->is           = NULL;
7760a545947SLisandro Dalcin     osm->is_local     = NULL;
77792bb6962SLisandro Dalcin     if (is) {
7789566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(n, &osm->is));
7792fa5cd67SKarl Rupp       for (i = 0; i < n; i++) osm->is[i] = is[i];
7803d03bb48SJed Brown       /* Flag indicating that the user has set overlapping subdomains so PCASM should not increase their size. */
7813d03bb48SJed Brown       osm->overlap = -1;
78292bb6962SLisandro Dalcin     }
7832b691e39SMatthew Knepley     if (is_local) {
7849566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(n, &osm->is_local));
7852fa5cd67SKarl Rupp       for (i = 0; i < n; i++) osm->is_local[i] = is_local[i];
786a35b7b57SDmitry Karpeev       if (!is) {
7879566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(osm->n_local_true, &osm->is));
788a35b7b57SDmitry Karpeev         for (i = 0; i < osm->n_local_true; i++) {
789a35b7b57SDmitry Karpeev           if (osm->overlap > 0) { /* With positive overlap, osm->is[i] will be modified */
7909566063dSJacob Faibussowitsch             PetscCall(ISDuplicate(osm->is_local[i], &osm->is[i]));
7919566063dSJacob Faibussowitsch             PetscCall(ISCopy(osm->is_local[i], osm->is[i]));
792a35b7b57SDmitry Karpeev           } else {
7939566063dSJacob Faibussowitsch             PetscCall(PetscObjectReference((PetscObject)osm->is_local[i]));
794a35b7b57SDmitry Karpeev             osm->is[i] = osm->is_local[i];
795a35b7b57SDmitry Karpeev           }
796a35b7b57SDmitry Karpeev         }
797a35b7b57SDmitry Karpeev       }
7982b691e39SMatthew Knepley     }
7994b9ad928SBarry Smith   }
8003ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8014b9ad928SBarry Smith }
8024b9ad928SBarry Smith 
803d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCASMSetTotalSubdomains_ASM(PC pc, PetscInt N, IS *is, IS *is_local)
804d71ae5a4SJacob Faibussowitsch {
8054b9ad928SBarry Smith   PC_ASM     *osm = (PC_ASM *)pc->data;
80613f74950SBarry Smith   PetscMPIInt rank, size;
80778904715SLisandro Dalcin   PetscInt    n;
8084b9ad928SBarry Smith 
8094b9ad928SBarry Smith   PetscFunctionBegin;
81063a3b9bcSJacob Faibussowitsch   PetscCheck(N >= 1, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Number of total blocks must be > 0, N = %" PetscInt_FMT, N);
81100045ab3SPierre Jolivet   PetscCheck(!is && !is_local, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Use PCASMSetLocalSubdomains() to set specific index sets, they cannot be set globally yet.");
8124b9ad928SBarry Smith 
8134b9ad928SBarry Smith   /*
814880770e9SJed Brown      Split the subdomains equally among all processors
8154b9ad928SBarry Smith   */
8169566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank));
8179566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size));
8184b9ad928SBarry Smith   n = N / size + ((N % size) > rank);
819835f2295SStefano Zampini   PetscCheck(n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Process %d must have at least one block: total processors %d total blocks %" PetscInt_FMT, rank, size, N);
8207827d75bSBarry Smith   PetscCheck(!pc->setupcalled || n == osm->n_local_true, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "PCASMSetTotalSubdomains() should be called before PCSetUp().");
8214b9ad928SBarry Smith   if (!pc->setupcalled) {
822ce78bad3SBarry Smith     PetscCall(PCASMDestroySubdomains(osm->n_local_true, &osm->is, &osm->is_local));
8232fa5cd67SKarl Rupp 
8244b9ad928SBarry Smith     osm->n_local_true = n;
8250a545947SLisandro Dalcin     osm->is           = NULL;
8260a545947SLisandro Dalcin     osm->is_local     = NULL;
8274b9ad928SBarry Smith   }
8283ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8294b9ad928SBarry Smith }
8304b9ad928SBarry Smith 
831d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCASMSetOverlap_ASM(PC pc, PetscInt ovl)
832d71ae5a4SJacob Faibussowitsch {
83392bb6962SLisandro Dalcin   PC_ASM *osm = (PC_ASM *)pc->data;
8344b9ad928SBarry Smith 
8354b9ad928SBarry Smith   PetscFunctionBegin;
8367827d75bSBarry Smith   PetscCheck(ovl >= 0, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Negative overlap value requested");
8377827d75bSBarry Smith   PetscCheck(!pc->setupcalled || ovl == osm->overlap, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "PCASMSetOverlap() should be called before PCSetUp().");
8382fa5cd67SKarl Rupp   if (!pc->setupcalled) osm->overlap = ovl;
8393ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8404b9ad928SBarry Smith }
8414b9ad928SBarry Smith 
842d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCASMSetType_ASM(PC pc, PCASMType type)
843d71ae5a4SJacob Faibussowitsch {
84492bb6962SLisandro Dalcin   PC_ASM *osm = (PC_ASM *)pc->data;
8454b9ad928SBarry Smith 
8464b9ad928SBarry Smith   PetscFunctionBegin;
8474b9ad928SBarry Smith   osm->type     = type;
848bf108f30SBarry Smith   osm->type_set = PETSC_TRUE;
8493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8504b9ad928SBarry Smith }
8514b9ad928SBarry Smith 
852d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCASMGetType_ASM(PC pc, PCASMType *type)
853d71ae5a4SJacob Faibussowitsch {
854c60c7ad4SBarry Smith   PC_ASM *osm = (PC_ASM *)pc->data;
855c60c7ad4SBarry Smith 
856c60c7ad4SBarry Smith   PetscFunctionBegin;
857c60c7ad4SBarry Smith   *type = osm->type;
8583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
859c60c7ad4SBarry Smith }
860c60c7ad4SBarry Smith 
861d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCASMSetLocalType_ASM(PC pc, PCCompositeType type)
862d71ae5a4SJacob Faibussowitsch {
86312cd4985SMatthew G. Knepley   PC_ASM *osm = (PC_ASM *)pc->data;
86412cd4985SMatthew G. Knepley 
86512cd4985SMatthew G. Knepley   PetscFunctionBegin;
8667827d75bSBarry 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");
86712cd4985SMatthew G. Knepley   osm->loctype = type;
8683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
86912cd4985SMatthew G. Knepley }
87012cd4985SMatthew G. Knepley 
871d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCASMGetLocalType_ASM(PC pc, PCCompositeType *type)
872d71ae5a4SJacob Faibussowitsch {
87312cd4985SMatthew G. Knepley   PC_ASM *osm = (PC_ASM *)pc->data;
87412cd4985SMatthew G. Knepley 
87512cd4985SMatthew G. Knepley   PetscFunctionBegin;
87612cd4985SMatthew G. Knepley   *type = osm->loctype;
8773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
87812cd4985SMatthew G. Knepley }
87912cd4985SMatthew G. Knepley 
880d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCASMSetSortIndices_ASM(PC pc, PetscBool doSort)
881d71ae5a4SJacob Faibussowitsch {
8826ed231c7SMatthew Knepley   PC_ASM *osm = (PC_ASM *)pc->data;
8836ed231c7SMatthew Knepley 
8846ed231c7SMatthew Knepley   PetscFunctionBegin;
8856ed231c7SMatthew Knepley   osm->sort_indices = doSort;
8863ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8876ed231c7SMatthew Knepley }
8886ed231c7SMatthew Knepley 
889d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCASMGetSubKSP_ASM(PC pc, PetscInt *n_local, PetscInt *first_local, KSP **ksp)
890d71ae5a4SJacob Faibussowitsch {
89192bb6962SLisandro Dalcin   PC_ASM *osm = (PC_ASM *)pc->data;
8924b9ad928SBarry Smith 
8934b9ad928SBarry Smith   PetscFunctionBegin;
8947827d75bSBarry 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");
8954b9ad928SBarry Smith 
8962fa5cd67SKarl Rupp   if (n_local) *n_local = osm->n_local_true;
89792bb6962SLisandro Dalcin   if (first_local) {
8989566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Scan(&osm->n_local_true, first_local, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)pc)));
89992bb6962SLisandro Dalcin     *first_local -= osm->n_local_true;
90092bb6962SLisandro Dalcin   }
901ed155784SPierre Jolivet   if (ksp) *ksp = osm->ksp;
9023ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9034b9ad928SBarry Smith }
9044b9ad928SBarry Smith 
905d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCASMGetSubMatType_ASM(PC pc, MatType *sub_mat_type)
906d71ae5a4SJacob Faibussowitsch {
90780ec0b7dSPatrick Sanan   PC_ASM *osm = (PC_ASM *)pc->data;
90880ec0b7dSPatrick Sanan 
90980ec0b7dSPatrick Sanan   PetscFunctionBegin;
91080ec0b7dSPatrick Sanan   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
9114f572ea9SToby Isaac   PetscAssertPointer(sub_mat_type, 2);
91280ec0b7dSPatrick Sanan   *sub_mat_type = osm->sub_mat_type;
9133ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
91480ec0b7dSPatrick Sanan }
91580ec0b7dSPatrick Sanan 
916d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCASMSetSubMatType_ASM(PC pc, MatType sub_mat_type)
917d71ae5a4SJacob Faibussowitsch {
91880ec0b7dSPatrick Sanan   PC_ASM *osm = (PC_ASM *)pc->data;
91980ec0b7dSPatrick Sanan 
92080ec0b7dSPatrick Sanan   PetscFunctionBegin;
92180ec0b7dSPatrick Sanan   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
9229566063dSJacob Faibussowitsch   PetscCall(PetscFree(osm->sub_mat_type));
9239566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(sub_mat_type, (char **)&osm->sub_mat_type));
9243ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
92580ec0b7dSPatrick Sanan }
92680ec0b7dSPatrick Sanan 
9275d83a8b1SBarry Smith /*@
928f1580f4eSBarry Smith   PCASMSetLocalSubdomains - Sets the local subdomains (for this processor only) for the additive Schwarz preconditioner `PCASM`.
9294b9ad928SBarry Smith 
930c3339decSBarry Smith   Collective
9314b9ad928SBarry Smith 
9324b9ad928SBarry Smith   Input Parameters:
9334b9ad928SBarry Smith + pc       - the preconditioner context
9344b9ad928SBarry Smith . n        - the number of subdomains for this processor (default value = 1)
935a3b724e8SBarry Smith . is       - the index set that defines the subdomains for this processor (or `NULL` for PETSc to determine subdomains)
936f13dfd9eSBarry Smith              the values of the `is` array are copied so you can free the array (not the `IS` in the array) after this call
937f13dfd9eSBarry 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`
938f13dfd9eSBarry Smith              (or `NULL` to not provide these). The values of the `is_local` array are copied so you can free the array
939f13dfd9eSBarry Smith              (not the `IS` in the array) after this call
9404b9ad928SBarry Smith 
941342c94f9SMatthew G. Knepley   Options Database Key:
94220f4b53cSBarry Smith . -pc_asm_local_blocks <blks> - Sets number of local blocks
94320f4b53cSBarry Smith 
94420f4b53cSBarry Smith   Level: advanced
945342c94f9SMatthew G. Knepley 
9464b9ad928SBarry Smith   Notes:
947a3b724e8SBarry Smith   The `IS` numbering is in the parallel, global numbering of the vector for both `is` and `is_local`
9484b9ad928SBarry Smith 
949f1580f4eSBarry Smith   By default the `PCASM` preconditioner uses 1 block per processor.
9504b9ad928SBarry Smith 
951f1580f4eSBarry Smith   Use `PCASMSetTotalSubdomains()` to set the subdomains for all processors.
9524b9ad928SBarry Smith 
953a3b724e8SBarry Smith   If `is_local` is provided and `PCASMType` is `PC_ASM_RESTRICT` then the solution only over the `is_local` region is interpolated
954a3b724e8SBarry Smith   back to form the global solution (this is the standard restricted additive Schwarz method, RASM)
955f1ee410cSBarry Smith 
956a3b724e8SBarry Smith   If `is_local` is provided and `PCASMType` is `PC_ASM_INTERPOLATE` or `PC_ASM_NONE` then an error is generated since there is
957f1ee410cSBarry Smith   no code to handle that case.
958f1ee410cSBarry Smith 
959562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMSetTotalSubdomains()`, `PCASMSetOverlap()`, `PCASMGetSubKSP()`,
960f1580f4eSBarry Smith           `PCASMCreateSubdomains2D()`, `PCASMGetLocalSubdomains()`, `PCASMType`, `PCASMSetType()`, `PCGASM`
9614b9ad928SBarry Smith @*/
962d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMSetLocalSubdomains(PC pc, PetscInt n, IS is[], IS is_local[])
963d71ae5a4SJacob Faibussowitsch {
9644b9ad928SBarry Smith   PetscFunctionBegin;
9650700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
966cac4c232SBarry Smith   PetscTryMethod(pc, "PCASMSetLocalSubdomains_C", (PC, PetscInt, IS[], IS[]), (pc, n, is, is_local));
9673ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9684b9ad928SBarry Smith }
9694b9ad928SBarry Smith 
9705d83a8b1SBarry Smith /*@
971feb221f8SDmitry Karpeev   PCASMSetTotalSubdomains - Sets the subdomains for all processors for the
972f1580f4eSBarry Smith   additive Schwarz preconditioner, `PCASM`.
9734b9ad928SBarry Smith 
974c3339decSBarry Smith   Collective, all MPI ranks must pass in the same array of `IS`
9754b9ad928SBarry Smith 
9764b9ad928SBarry Smith   Input Parameters:
9774b9ad928SBarry Smith + pc       - the preconditioner context
978feb221f8SDmitry Karpeev . N        - the number of subdomains for all processors
979a3b724e8SBarry Smith . is       - the index sets that define the subdomains for all processors (or `NULL` to ask PETSc to determine the subdomains)
9805d83a8b1SBarry Smith              the values of the `is` array are copied so you can free the array (not the `IS` in the array) after this call
981a3b724e8SBarry Smith - is_local - the index sets that define the local part of the subdomains for this processor (or `NULL` to not provide this information)
9825d83a8b1SBarry Smith              The values of the `is_local` array are copied so you can free the array (not the `IS` in the array) after this call
9834b9ad928SBarry Smith 
9844b9ad928SBarry Smith   Options Database Key:
9854b9ad928SBarry Smith . -pc_asm_blocks <blks> - Sets total blocks
9864b9ad928SBarry Smith 
98720f4b53cSBarry Smith   Level: advanced
98820f4b53cSBarry Smith 
9894b9ad928SBarry Smith   Notes:
990a3b724e8SBarry Smith   Currently you cannot use this to set the actual subdomains with the argument `is` or `is_local`.
9914b9ad928SBarry Smith 
992f1580f4eSBarry Smith   By default the `PCASM` preconditioner uses 1 block per processor.
9934b9ad928SBarry Smith 
9944b9ad928SBarry Smith   These index sets cannot be destroyed until after completion of the
995f1580f4eSBarry Smith   linear solves for which the `PCASM` preconditioner is being used.
9964b9ad928SBarry Smith 
997f1580f4eSBarry Smith   Use `PCASMSetLocalSubdomains()` to set local subdomains.
9984b9ad928SBarry Smith 
999f1580f4eSBarry Smith   The `IS` numbering is in the parallel, global numbering of the vector for both is and is_local
10001093a601SBarry Smith 
1001562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMSetLocalSubdomains()`, `PCASMSetOverlap()`, `PCASMGetSubKSP()`,
1002f1580f4eSBarry Smith           `PCASMCreateSubdomains2D()`, `PCGASM`
10034b9ad928SBarry Smith @*/
1004d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMSetTotalSubdomains(PC pc, PetscInt N, IS is[], IS is_local[])
1005d71ae5a4SJacob Faibussowitsch {
10064b9ad928SBarry Smith   PetscFunctionBegin;
10070700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1008cac4c232SBarry Smith   PetscTryMethod(pc, "PCASMSetTotalSubdomains_C", (PC, PetscInt, IS[], IS[]), (pc, N, is, is_local));
10093ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10104b9ad928SBarry Smith }
10114b9ad928SBarry Smith 
10124b9ad928SBarry Smith /*@
10134b9ad928SBarry Smith   PCASMSetOverlap - Sets the overlap between a pair of subdomains for the
1014f1580f4eSBarry Smith   additive Schwarz preconditioner, `PCASM`.
10154b9ad928SBarry Smith 
1016c3339decSBarry Smith   Logically Collective
10174b9ad928SBarry Smith 
10184b9ad928SBarry Smith   Input Parameters:
10194b9ad928SBarry Smith + pc  - the preconditioner context
10204b9ad928SBarry Smith - ovl - the amount of overlap between subdomains (ovl >= 0, default value = 1)
10214b9ad928SBarry Smith 
10224b9ad928SBarry Smith   Options Database Key:
10234b9ad928SBarry Smith . -pc_asm_overlap <ovl> - Sets overlap
10244b9ad928SBarry Smith 
102520f4b53cSBarry Smith   Level: intermediate
102620f4b53cSBarry Smith 
10274b9ad928SBarry Smith   Notes:
1028f1580f4eSBarry Smith   By default the `PCASM` preconditioner uses 1 block per processor.  To use
1029f1580f4eSBarry Smith   multiple blocks per perocessor, see `PCASMSetTotalSubdomains()` and
1030f1580f4eSBarry Smith   `PCASMSetLocalSubdomains()` (and the option -pc_asm_blocks <blks>).
10314b9ad928SBarry Smith 
10324b9ad928SBarry Smith   The overlap defaults to 1, so if one desires that no additional
10334b9ad928SBarry Smith   overlap be computed beyond what may have been set with a call to
1034f1580f4eSBarry Smith   `PCASMSetTotalSubdomains()` or `PCASMSetLocalSubdomains()`, then ovl
10354b9ad928SBarry Smith   must be set to be 0.  In particular, if one does not explicitly set
10364b9ad928SBarry Smith   the subdomains an application code, then all overlap would be computed
1037f1580f4eSBarry Smith   internally by PETSc, and using an overlap of 0 would result in an `PCASM`
10384b9ad928SBarry Smith   variant that is equivalent to the block Jacobi preconditioner.
10394b9ad928SBarry Smith 
1040f1ee410cSBarry Smith   The default algorithm used by PETSc to increase overlap is fast, but not scalable,
1041f1ee410cSBarry Smith   use the option -mat_increase_overlap_scalable when the problem and number of processes is large.
1042f1ee410cSBarry Smith 
10432fe279fdSBarry Smith   One can define initial index sets with any overlap via
1044f1580f4eSBarry Smith   `PCASMSetLocalSubdomains()`; the routine
1045f1580f4eSBarry Smith   `PCASMSetOverlap()` merely allows PETSc to extend that overlap further
10464b9ad928SBarry Smith   if desired.
10474b9ad928SBarry Smith 
1048562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMSetTotalSubdomains()`, `PCASMSetLocalSubdomains()`, `PCASMGetSubKSP()`,
1049f1580f4eSBarry Smith           `PCASMCreateSubdomains2D()`, `PCASMGetLocalSubdomains()`, `MatIncreaseOverlap()`, `PCGASM`
10504b9ad928SBarry Smith @*/
1051d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMSetOverlap(PC pc, PetscInt ovl)
1052d71ae5a4SJacob Faibussowitsch {
10534b9ad928SBarry Smith   PetscFunctionBegin;
10540700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1055c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc, ovl, 2);
1056cac4c232SBarry Smith   PetscTryMethod(pc, "PCASMSetOverlap_C", (PC, PetscInt), (pc, ovl));
10573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10584b9ad928SBarry Smith }
10594b9ad928SBarry Smith 
10604b9ad928SBarry Smith /*@
10614b9ad928SBarry Smith   PCASMSetType - Sets the type of restriction and interpolation used
1062f1580f4eSBarry Smith   for local problems in the additive Schwarz method, `PCASM`.
10634b9ad928SBarry Smith 
1064c3339decSBarry Smith   Logically Collective
10654b9ad928SBarry Smith 
10664b9ad928SBarry Smith   Input Parameters:
10674b9ad928SBarry Smith + pc   - the preconditioner context
1068f1580f4eSBarry Smith - type - variant of `PCASM`, one of
10694b9ad928SBarry Smith .vb
10704b9ad928SBarry Smith       PC_ASM_BASIC       - full interpolation and restriction
107182b5ce2aSStefano Zampini       PC_ASM_RESTRICT    - full restriction, local processor interpolation (default)
10724b9ad928SBarry Smith       PC_ASM_INTERPOLATE - full interpolation, local processor restriction
10734b9ad928SBarry Smith       PC_ASM_NONE        - local processor restriction and interpolation
10744b9ad928SBarry Smith .ve
10754b9ad928SBarry Smith 
10764b9ad928SBarry Smith   Options Database Key:
1077f1580f4eSBarry Smith . -pc_asm_type [basic,restrict,interpolate,none] - Sets `PCASMType`
10784b9ad928SBarry Smith 
107920f4b53cSBarry Smith   Level: intermediate
108020f4b53cSBarry Smith 
1081f1580f4eSBarry Smith   Note:
1082f1580f4eSBarry Smith   if the is_local arguments are passed to `PCASMSetLocalSubdomains()` then they are used when `PC_ASM_RESTRICT` has been selected
1083f1ee410cSBarry Smith   to limit the local processor interpolation
1084f1ee410cSBarry Smith 
1085562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMSetTotalSubdomains()`, `PCASMGetSubKSP()`,
1086f1580f4eSBarry Smith           `PCASMCreateSubdomains2D()`, `PCASMType`, `PCASMSetLocalType()`, `PCASMGetLocalType()`, `PCGASM`
10874b9ad928SBarry Smith @*/
1088d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMSetType(PC pc, PCASMType type)
1089d71ae5a4SJacob Faibussowitsch {
10904b9ad928SBarry Smith   PetscFunctionBegin;
10910700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1092c5eb9154SBarry Smith   PetscValidLogicalCollectiveEnum(pc, type, 2);
1093cac4c232SBarry Smith   PetscTryMethod(pc, "PCASMSetType_C", (PC, PCASMType), (pc, type));
10943ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10954b9ad928SBarry Smith }
10964b9ad928SBarry Smith 
1097c60c7ad4SBarry Smith /*@
1098c60c7ad4SBarry Smith   PCASMGetType - Gets the type of restriction and interpolation used
1099f1580f4eSBarry Smith   for local problems in the additive Schwarz method, `PCASM`.
1100c60c7ad4SBarry Smith 
1101c3339decSBarry Smith   Logically Collective
1102c60c7ad4SBarry Smith 
1103c60c7ad4SBarry Smith   Input Parameter:
1104c60c7ad4SBarry Smith . pc - the preconditioner context
1105c60c7ad4SBarry Smith 
1106c60c7ad4SBarry Smith   Output Parameter:
1107f1580f4eSBarry Smith . type - variant of `PCASM`, one of
1108c60c7ad4SBarry Smith .vb
1109c60c7ad4SBarry Smith       PC_ASM_BASIC       - full interpolation and restriction
1110c60c7ad4SBarry Smith       PC_ASM_RESTRICT    - full restriction, local processor interpolation
1111c60c7ad4SBarry Smith       PC_ASM_INTERPOLATE - full interpolation, local processor restriction
1112c60c7ad4SBarry Smith       PC_ASM_NONE        - local processor restriction and interpolation
1113c60c7ad4SBarry Smith .ve
1114c60c7ad4SBarry Smith 
1115c60c7ad4SBarry Smith   Options Database Key:
1116f1580f4eSBarry Smith . -pc_asm_type [basic,restrict,interpolate,none] - Sets `PCASM` type
1117c60c7ad4SBarry Smith 
1118c60c7ad4SBarry Smith   Level: intermediate
1119c60c7ad4SBarry Smith 
1120562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMSetTotalSubdomains()`, `PCASMGetSubKSP()`, `PCGASM`,
1121db781477SPatrick Sanan           `PCASMCreateSubdomains2D()`, `PCASMType`, `PCASMSetType()`, `PCASMSetLocalType()`, `PCASMGetLocalType()`
1122c60c7ad4SBarry Smith @*/
1123d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMGetType(PC pc, PCASMType *type)
1124d71ae5a4SJacob Faibussowitsch {
1125c60c7ad4SBarry Smith   PetscFunctionBegin;
1126c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1127cac4c232SBarry Smith   PetscUseMethod(pc, "PCASMGetType_C", (PC, PCASMType *), (pc, type));
11283ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1129c60c7ad4SBarry Smith }
1130c60c7ad4SBarry Smith 
113112cd4985SMatthew G. Knepley /*@
1132f1580f4eSBarry Smith   PCASMSetLocalType - Sets the type of composition used for local problems in the additive Schwarz method, `PCASM`.
113312cd4985SMatthew G. Knepley 
1134c3339decSBarry Smith   Logically Collective
113512cd4985SMatthew G. Knepley 
113612cd4985SMatthew G. Knepley   Input Parameters:
113712cd4985SMatthew G. Knepley + pc   - the preconditioner context
113812cd4985SMatthew G. Knepley - type - type of composition, one of
113912cd4985SMatthew G. Knepley .vb
114012cd4985SMatthew G. Knepley   PC_COMPOSITE_ADDITIVE       - local additive combination
114112cd4985SMatthew G. Knepley   PC_COMPOSITE_MULTIPLICATIVE - local multiplicative combination
114212cd4985SMatthew G. Knepley .ve
114312cd4985SMatthew G. Knepley 
114412cd4985SMatthew G. Knepley   Options Database Key:
114512cd4985SMatthew G. Knepley . -pc_asm_local_type [additive,multiplicative] - Sets local solver composition type
114612cd4985SMatthew G. Knepley 
114712cd4985SMatthew G. Knepley   Level: intermediate
114812cd4985SMatthew G. Knepley 
1149562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMSetType()`, `PCASMGetType()`, `PCASMGetLocalType()`, `PCASMType`, `PCCompositeType`
115012cd4985SMatthew G. Knepley @*/
1151d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMSetLocalType(PC pc, PCCompositeType type)
1152d71ae5a4SJacob Faibussowitsch {
115312cd4985SMatthew G. Knepley   PetscFunctionBegin;
115412cd4985SMatthew G. Knepley   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
115512cd4985SMatthew G. Knepley   PetscValidLogicalCollectiveEnum(pc, type, 2);
1156cac4c232SBarry Smith   PetscTryMethod(pc, "PCASMSetLocalType_C", (PC, PCCompositeType), (pc, type));
11573ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
115812cd4985SMatthew G. Knepley }
115912cd4985SMatthew G. Knepley 
116012cd4985SMatthew G. Knepley /*@
1161f1580f4eSBarry Smith   PCASMGetLocalType - Gets the type of composition used for local problems in the additive Schwarz method, `PCASM`.
116212cd4985SMatthew G. Knepley 
1163c3339decSBarry Smith   Logically Collective
116412cd4985SMatthew G. Knepley 
116512cd4985SMatthew G. Knepley   Input Parameter:
116612cd4985SMatthew G. Knepley . pc - the preconditioner context
116712cd4985SMatthew G. Knepley 
116812cd4985SMatthew G. Knepley   Output Parameter:
116912cd4985SMatthew G. Knepley . type - type of composition, one of
117012cd4985SMatthew G. Knepley .vb
117112cd4985SMatthew G. Knepley   PC_COMPOSITE_ADDITIVE       - local additive combination
117212cd4985SMatthew G. Knepley   PC_COMPOSITE_MULTIPLICATIVE - local multiplicative combination
117312cd4985SMatthew G. Knepley .ve
117412cd4985SMatthew G. Knepley 
117512cd4985SMatthew G. Knepley   Options Database Key:
117612cd4985SMatthew G. Knepley . -pc_asm_local_type [additive,multiplicative] - Sets local solver composition type
117712cd4985SMatthew G. Knepley 
117812cd4985SMatthew G. Knepley   Level: intermediate
117912cd4985SMatthew G. Knepley 
11800241b274SPierre Jolivet .seealso: [](ch_ksp), `PCASM`, `PCASMSetType()`, `PCASMGetType()`, `PCASMSetLocalType()`, `PCASMType`, `PCCompositeType`
118112cd4985SMatthew G. Knepley @*/
1182d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMGetLocalType(PC pc, PCCompositeType *type)
1183d71ae5a4SJacob Faibussowitsch {
118412cd4985SMatthew G. Knepley   PetscFunctionBegin;
118512cd4985SMatthew G. Knepley   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
11864f572ea9SToby Isaac   PetscAssertPointer(type, 2);
1187cac4c232SBarry Smith   PetscUseMethod(pc, "PCASMGetLocalType_C", (PC, PCCompositeType *), (pc, type));
11883ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
118912cd4985SMatthew G. Knepley }
119012cd4985SMatthew G. Knepley 
11916ed231c7SMatthew Knepley /*@
11926ed231c7SMatthew Knepley   PCASMSetSortIndices - Determines whether subdomain indices are sorted.
11936ed231c7SMatthew Knepley 
1194c3339decSBarry Smith   Logically Collective
11956ed231c7SMatthew Knepley 
11966ed231c7SMatthew Knepley   Input Parameters:
11976ed231c7SMatthew Knepley + pc     - the preconditioner context
11986ed231c7SMatthew Knepley - doSort - sort the subdomain indices
11996ed231c7SMatthew Knepley 
12006ed231c7SMatthew Knepley   Level: intermediate
12016ed231c7SMatthew Knepley 
1202562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMSetLocalSubdomains()`, `PCASMSetTotalSubdomains()`, `PCASMGetSubKSP()`,
1203db781477SPatrick Sanan           `PCASMCreateSubdomains2D()`
12046ed231c7SMatthew Knepley @*/
1205d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMSetSortIndices(PC pc, PetscBool doSort)
1206d71ae5a4SJacob Faibussowitsch {
12076ed231c7SMatthew Knepley   PetscFunctionBegin;
12080700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1209acfcf0e5SJed Brown   PetscValidLogicalCollectiveBool(pc, doSort, 2);
1210cac4c232SBarry Smith   PetscTryMethod(pc, "PCASMSetSortIndices_C", (PC, PetscBool), (pc, doSort));
12113ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
12126ed231c7SMatthew Knepley }
12136ed231c7SMatthew Knepley 
12144b9ad928SBarry Smith /*@C
1215f1580f4eSBarry Smith   PCASMGetSubKSP - Gets the local `KSP` contexts for all blocks on
12164b9ad928SBarry Smith   this processor.
12174b9ad928SBarry Smith 
1218c3339decSBarry Smith   Collective iff first_local is requested
12194b9ad928SBarry Smith 
12204b9ad928SBarry Smith   Input Parameter:
12214b9ad928SBarry Smith . pc - the preconditioner context
12224b9ad928SBarry Smith 
12234b9ad928SBarry Smith   Output Parameters:
1224a3b724e8SBarry Smith + n_local     - the number of blocks on this processor or `NULL`
1225a3b724e8SBarry Smith . first_local - the global number of the first block on this processor or `NULL`, all processors must request or all must pass `NULL`
1226f1580f4eSBarry Smith - ksp         - the array of `KSP` contexts
12274b9ad928SBarry Smith 
122820f4b53cSBarry Smith   Level: advanced
122920f4b53cSBarry Smith 
1230f1580f4eSBarry Smith   Notes:
1231f1580f4eSBarry Smith   After `PCASMGetSubKSP()` the array of `KSP`s is not to be freed.
12324b9ad928SBarry Smith 
1233f1580f4eSBarry Smith   You must call `KSPSetUp()` before calling `PCASMGetSubKSP()`.
12344b9ad928SBarry Smith 
1235562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMSetTotalSubdomains()`, `PCASMSetOverlap()`,
1236db781477SPatrick Sanan           `PCASMCreateSubdomains2D()`,
12374b9ad928SBarry Smith @*/
1238d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMGetSubKSP(PC pc, PetscInt *n_local, PetscInt *first_local, KSP *ksp[])
1239d71ae5a4SJacob Faibussowitsch {
12404b9ad928SBarry Smith   PetscFunctionBegin;
12410700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1242cac4c232SBarry Smith   PetscUseMethod(pc, "PCASMGetSubKSP_C", (PC, PetscInt *, PetscInt *, KSP **), (pc, n_local, first_local, ksp));
12433ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
12444b9ad928SBarry Smith }
12454b9ad928SBarry Smith 
12464b9ad928SBarry Smith /*MC
12473b09bd56SBarry Smith    PCASM - Use the (restricted) additive Schwarz method, each block is (approximately) solved with
12481d27aa22SBarry Smith            its own `KSP` object, {cite}`dryja1987additive` and {cite}`1sbg`
12494b9ad928SBarry Smith 
12504b9ad928SBarry Smith    Options Database Keys:
12511d27aa22SBarry Smith +  -pc_asm_blocks <blks>                          - Sets total blocks. Defaults to one block per MPI process.
12524b9ad928SBarry Smith .  -pc_asm_overlap <ovl>                          - Sets overlap
1253f1580f4eSBarry Smith .  -pc_asm_type [basic,restrict,interpolate,none] - Sets `PCASMType`, default is restrict. See `PCASMSetType()`
125473ff1848SBarry Smith .  -pc_asm_dm_subdomains <bool>                   - use subdomains defined by the `DM` with `DMCreateDomainDecomposition()`
1255f1580f4eSBarry Smith -  -pc_asm_local_type [additive, multiplicative]  - Sets `PCCompositeType`, default is additive. See `PCASMSetLocalType()`
12564b9ad928SBarry Smith 
12574b9ad928SBarry Smith    Level: beginner
12584b9ad928SBarry Smith 
1259f1580f4eSBarry Smith    Notes:
1260f1580f4eSBarry Smith    If you run with, for example, 3 blocks on 1 processor or 3 blocks on 3 processors you
126173ff1848SBarry Smith    will get a different convergence rate due to the default option of `-pc_asm_type restrict`. Use
126273ff1848SBarry Smith    `-pc_asm_type basic` to get the same convergence behavior
1263f1580f4eSBarry Smith 
1264f1580f4eSBarry Smith    Each processor can have one or more blocks, but a block cannot be shared by more
1265f1580f4eSBarry Smith    than one processor. Use `PCGASM` for subdomains shared by multiple processes.
1266f1580f4eSBarry Smith 
126773ff1848SBarry Smith    To set options on the solvers for each block append `-sub_` to all the `KSP`, and `PC`
126873ff1848SBarry Smith    options database keys. For example, `-sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly`
1269f1580f4eSBarry Smith 
1270f1580f4eSBarry Smith    To set the options on the solvers separate for each block call `PCASMGetSubKSP()`
1271f1580f4eSBarry Smith    and set the options directly on the resulting `KSP` object (you can access its `PC` with `KSPGetPC()`)
1272f1580f4eSBarry Smith 
127373ff1848SBarry Smith    If the `PC` has an associated `DM`, then, by default, `DMCreateDomainDecomposition()` is used to create the subdomains
127473ff1848SBarry Smith 
1275562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCASMType`, `PCCompositeType`,
1276db781477SPatrick Sanan           `PCBJACOBI`, `PCASMGetSubKSP()`, `PCASMSetLocalSubdomains()`, `PCASMType`, `PCASMGetType()`, `PCASMSetLocalType()`, `PCASMGetLocalType()`
1277db781477SPatrick Sanan           `PCASMSetTotalSubdomains()`, `PCSetModifySubMatrices()`, `PCASMSetOverlap()`, `PCASMSetType()`, `PCCompositeType`
12784b9ad928SBarry Smith M*/
12794b9ad928SBarry Smith 
1280d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_ASM(PC pc)
1281d71ae5a4SJacob Faibussowitsch {
12824b9ad928SBarry Smith   PC_ASM *osm;
12834b9ad928SBarry Smith 
12844b9ad928SBarry Smith   PetscFunctionBegin;
12854dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&osm));
12862fa5cd67SKarl Rupp 
12874b9ad928SBarry Smith   osm->n             = PETSC_DECIDE;
12884b9ad928SBarry Smith   osm->n_local       = 0;
12892b837212SDmitry Karpeev   osm->n_local_true  = PETSC_DECIDE;
12904b9ad928SBarry Smith   osm->overlap       = 1;
12910a545947SLisandro Dalcin   osm->ksp           = NULL;
12920a545947SLisandro Dalcin   osm->restriction   = NULL;
12930a545947SLisandro Dalcin   osm->lprolongation = NULL;
12940a545947SLisandro Dalcin   osm->lrestriction  = NULL;
12950a545947SLisandro Dalcin   osm->x             = NULL;
12960a545947SLisandro Dalcin   osm->y             = NULL;
12970a545947SLisandro Dalcin   osm->is            = NULL;
12980a545947SLisandro Dalcin   osm->is_local      = NULL;
12990a545947SLisandro Dalcin   osm->mat           = NULL;
13000a545947SLisandro Dalcin   osm->pmat          = NULL;
13014b9ad928SBarry Smith   osm->type          = PC_ASM_RESTRICT;
130212cd4985SMatthew G. Knepley   osm->loctype       = PC_COMPOSITE_ADDITIVE;
13036ed231c7SMatthew Knepley   osm->sort_indices  = PETSC_TRUE;
1304d709ea83SDmitry Karpeev   osm->dm_subdomains = PETSC_FALSE;
130580ec0b7dSPatrick Sanan   osm->sub_mat_type  = NULL;
13064b9ad928SBarry Smith 
130792bb6962SLisandro Dalcin   pc->data                   = (void *)osm;
13084b9ad928SBarry Smith   pc->ops->apply             = PCApply_ASM;
130948e38a8aSPierre Jolivet   pc->ops->matapply          = PCMatApply_ASM;
13104b9ad928SBarry Smith   pc->ops->applytranspose    = PCApplyTranspose_ASM;
1311*4dbf25a8SPierre Jolivet   pc->ops->matapplytranspose = PCMatApplyTranspose_ASM;
13124b9ad928SBarry Smith   pc->ops->setup             = PCSetUp_ASM;
1313e91c6855SBarry Smith   pc->ops->reset             = PCReset_ASM;
13144b9ad928SBarry Smith   pc->ops->destroy           = PCDestroy_ASM;
13154b9ad928SBarry Smith   pc->ops->setfromoptions    = PCSetFromOptions_ASM;
13164b9ad928SBarry Smith   pc->ops->setuponblocks     = PCSetUpOnBlocks_ASM;
13174b9ad928SBarry Smith   pc->ops->view              = PCView_ASM;
13180a545947SLisandro Dalcin   pc->ops->applyrichardson   = NULL;
13194b9ad928SBarry Smith 
13209566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetLocalSubdomains_C", PCASMSetLocalSubdomains_ASM));
13219566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetTotalSubdomains_C", PCASMSetTotalSubdomains_ASM));
13229566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetOverlap_C", PCASMSetOverlap_ASM));
13239566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetType_C", PCASMSetType_ASM));
13249566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMGetType_C", PCASMGetType_ASM));
13259566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetLocalType_C", PCASMSetLocalType_ASM));
13269566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMGetLocalType_C", PCASMGetLocalType_ASM));
13279566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetSortIndices_C", PCASMSetSortIndices_ASM));
13289566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMGetSubKSP_C", PCASMGetSubKSP_ASM));
13299566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMGetSubMatType_C", PCASMGetSubMatType_ASM));
13309566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetSubMatType_C", PCASMSetSubMatType_ASM));
13313ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13324b9ad928SBarry Smith }
13334b9ad928SBarry Smith 
133492bb6962SLisandro Dalcin /*@C
133592bb6962SLisandro Dalcin   PCASMCreateSubdomains - Creates the index sets for the overlapping Schwarz
1336f1580f4eSBarry Smith   preconditioner, `PCASM`,  for any problem on a general grid.
133792bb6962SLisandro Dalcin 
133892bb6962SLisandro Dalcin   Collective
133992bb6962SLisandro Dalcin 
134092bb6962SLisandro Dalcin   Input Parameters:
134192bb6962SLisandro Dalcin + A - The global matrix operator
134292bb6962SLisandro Dalcin - n - the number of local blocks
134392bb6962SLisandro Dalcin 
13442fe279fdSBarry Smith   Output Parameter:
134592bb6962SLisandro Dalcin . outis - the array of index sets defining the subdomains
134692bb6962SLisandro Dalcin 
134792bb6962SLisandro Dalcin   Level: advanced
134892bb6962SLisandro Dalcin 
1349f1580f4eSBarry Smith   Note:
1350f1580f4eSBarry Smith   This generates nonoverlapping subdomains; the `PCASM` will generate the overlap
1351f1580f4eSBarry Smith   from these if you use `PCASMSetLocalSubdomains()`
13527d6bfa3bSBarry Smith 
1353562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMSetLocalSubdomains()`, `PCASMDestroySubdomains()`
135492bb6962SLisandro Dalcin @*/
1355d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMCreateSubdomains(Mat A, PetscInt n, IS *outis[])
1356d71ae5a4SJacob Faibussowitsch {
135792bb6962SLisandro Dalcin   MatPartitioning mpart;
135892bb6962SLisandro Dalcin   const char     *prefix;
135992bb6962SLisandro Dalcin   PetscInt        i, j, rstart, rend, bs;
1360976e8c5aSLisandro Dalcin   PetscBool       hasop, isbaij = PETSC_FALSE, foundpart = PETSC_FALSE;
13610298fd71SBarry Smith   Mat             Ad = NULL, adj;
136292bb6962SLisandro Dalcin   IS              ispart, isnumb, *is;
136392bb6962SLisandro Dalcin 
136492bb6962SLisandro Dalcin   PetscFunctionBegin;
13650700a824SBarry Smith   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
13664f572ea9SToby Isaac   PetscAssertPointer(outis, 3);
136763a3b9bcSJacob Faibussowitsch   PetscCheck(n >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "number of local blocks must be > 0, n = %" PetscInt_FMT, n);
136892bb6962SLisandro Dalcin 
136992bb6962SLisandro Dalcin   /* Get prefix, row distribution, and block size */
13709566063dSJacob Faibussowitsch   PetscCall(MatGetOptionsPrefix(A, &prefix));
13719566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(A, &rstart, &rend));
13729566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(A, &bs));
137363a3b9bcSJacob 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);
137465e19b50SBarry Smith 
137592bb6962SLisandro Dalcin   /* Get diagonal block from matrix if possible */
13769566063dSJacob Faibussowitsch   PetscCall(MatHasOperation(A, MATOP_GET_DIAGONAL_BLOCK, &hasop));
137748a46eb9SPierre Jolivet   if (hasop) PetscCall(MatGetDiagonalBlock(A, &Ad));
137892bb6962SLisandro Dalcin   if (Ad) {
13799566063dSJacob Faibussowitsch     PetscCall(PetscObjectBaseTypeCompare((PetscObject)Ad, MATSEQBAIJ, &isbaij));
13809566063dSJacob Faibussowitsch     if (!isbaij) PetscCall(PetscObjectBaseTypeCompare((PetscObject)Ad, MATSEQSBAIJ, &isbaij));
138192bb6962SLisandro Dalcin   }
138292bb6962SLisandro Dalcin   if (Ad && n > 1) {
1383ace3abfcSBarry Smith     PetscBool match, done;
138492bb6962SLisandro Dalcin     /* Try to setup a good matrix partitioning if available */
13859566063dSJacob Faibussowitsch     PetscCall(MatPartitioningCreate(PETSC_COMM_SELF, &mpart));
13869566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)mpart, prefix));
13879566063dSJacob Faibussowitsch     PetscCall(MatPartitioningSetFromOptions(mpart));
13889566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)mpart, MATPARTITIONINGCURRENT, &match));
138948a46eb9SPierre Jolivet     if (!match) PetscCall(PetscObjectTypeCompare((PetscObject)mpart, MATPARTITIONINGSQUARE, &match));
139092bb6962SLisandro Dalcin     if (!match) { /* assume a "good" partitioner is available */
13911a83f524SJed Brown       PetscInt        na;
13921a83f524SJed Brown       const PetscInt *ia, *ja;
13939566063dSJacob Faibussowitsch       PetscCall(MatGetRowIJ(Ad, 0, PETSC_TRUE, isbaij, &na, &ia, &ja, &done));
139492bb6962SLisandro Dalcin       if (done) {
139592bb6962SLisandro Dalcin         /* Build adjacency matrix by hand. Unfortunately a call to
139692bb6962SLisandro Dalcin            MatConvert(Ad,MATMPIADJ,MAT_INITIAL_MATRIX,&adj) will
139792bb6962SLisandro Dalcin            remove the block-aij structure and we cannot expect
139892bb6962SLisandro Dalcin            MatPartitioning to split vertices as we need */
13990a545947SLisandro Dalcin         PetscInt        i, j, len, nnz, cnt, *iia = NULL, *jja = NULL;
14001a83f524SJed Brown         const PetscInt *row;
140192bb6962SLisandro Dalcin         nnz = 0;
140292bb6962SLisandro Dalcin         for (i = 0; i < na; i++) { /* count number of nonzeros */
140392bb6962SLisandro Dalcin           len = ia[i + 1] - ia[i];
140492bb6962SLisandro Dalcin           row = ja + ia[i];
140592bb6962SLisandro Dalcin           for (j = 0; j < len; j++) {
140692bb6962SLisandro Dalcin             if (row[j] == i) { /* don't count diagonal */
14079371c9d4SSatish Balay               len--;
14089371c9d4SSatish Balay               break;
140992bb6962SLisandro Dalcin             }
141092bb6962SLisandro Dalcin           }
141192bb6962SLisandro Dalcin           nnz += len;
141292bb6962SLisandro Dalcin         }
14139566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(na + 1, &iia));
14149566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(nnz, &jja));
141592bb6962SLisandro Dalcin         nnz    = 0;
141692bb6962SLisandro Dalcin         iia[0] = 0;
141792bb6962SLisandro Dalcin         for (i = 0; i < na; i++) { /* fill adjacency */
141892bb6962SLisandro Dalcin           cnt = 0;
141992bb6962SLisandro Dalcin           len = ia[i + 1] - ia[i];
142092bb6962SLisandro Dalcin           row = ja + ia[i];
142192bb6962SLisandro Dalcin           for (j = 0; j < len; j++) {
142292bb6962SLisandro Dalcin             if (row[j] != i) { /* if not diagonal */
142392bb6962SLisandro Dalcin               jja[nnz + cnt++] = row[j];
142492bb6962SLisandro Dalcin             }
142592bb6962SLisandro Dalcin           }
142692bb6962SLisandro Dalcin           nnz += cnt;
142792bb6962SLisandro Dalcin           iia[i + 1] = nnz;
142892bb6962SLisandro Dalcin         }
142992bb6962SLisandro Dalcin         /* Partitioning of the adjacency matrix */
14309566063dSJacob Faibussowitsch         PetscCall(MatCreateMPIAdj(PETSC_COMM_SELF, na, na, iia, jja, NULL, &adj));
14319566063dSJacob Faibussowitsch         PetscCall(MatPartitioningSetAdjacency(mpart, adj));
14329566063dSJacob Faibussowitsch         PetscCall(MatPartitioningSetNParts(mpart, n));
14339566063dSJacob Faibussowitsch         PetscCall(MatPartitioningApply(mpart, &ispart));
14349566063dSJacob Faibussowitsch         PetscCall(ISPartitioningToNumbering(ispart, &isnumb));
14359566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&adj));
143692bb6962SLisandro Dalcin         foundpart = PETSC_TRUE;
143792bb6962SLisandro Dalcin       }
14389566063dSJacob Faibussowitsch       PetscCall(MatRestoreRowIJ(Ad, 0, PETSC_TRUE, isbaij, &na, &ia, &ja, &done));
143992bb6962SLisandro Dalcin     }
14409566063dSJacob Faibussowitsch     PetscCall(MatPartitioningDestroy(&mpart));
144192bb6962SLisandro Dalcin   }
144292bb6962SLisandro Dalcin 
14439566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n, &is));
144492bb6962SLisandro Dalcin   *outis = is;
144592bb6962SLisandro Dalcin 
144692bb6962SLisandro Dalcin   if (!foundpart) {
144792bb6962SLisandro Dalcin     /* Partitioning by contiguous chunks of rows */
144892bb6962SLisandro Dalcin 
144992bb6962SLisandro Dalcin     PetscInt mbs   = (rend - rstart) / bs;
145092bb6962SLisandro Dalcin     PetscInt start = rstart;
145192bb6962SLisandro Dalcin     for (i = 0; i < n; i++) {
145292bb6962SLisandro Dalcin       PetscInt count = (mbs / n + ((mbs % n) > i)) * bs;
14539566063dSJacob Faibussowitsch       PetscCall(ISCreateStride(PETSC_COMM_SELF, count, start, 1, &is[i]));
145492bb6962SLisandro Dalcin       start += count;
145592bb6962SLisandro Dalcin     }
145692bb6962SLisandro Dalcin 
145792bb6962SLisandro Dalcin   } else {
145892bb6962SLisandro Dalcin     /* Partitioning by adjacency of diagonal block  */
145992bb6962SLisandro Dalcin 
146092bb6962SLisandro Dalcin     const PetscInt *numbering;
146192bb6962SLisandro Dalcin     PetscInt       *count, nidx, *indices, *newidx, start = 0;
146292bb6962SLisandro Dalcin     /* Get node count in each partition */
14639566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &count));
14649566063dSJacob Faibussowitsch     PetscCall(ISPartitioningCount(ispart, n, count));
146592bb6962SLisandro Dalcin     if (isbaij && bs > 1) { /* adjust for the block-aij case */
146692bb6962SLisandro Dalcin       for (i = 0; i < n; i++) count[i] *= bs;
146792bb6962SLisandro Dalcin     }
146892bb6962SLisandro Dalcin     /* Build indices from node numbering */
14699566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(isnumb, &nidx));
14709566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nidx, &indices));
147192bb6962SLisandro Dalcin     for (i = 0; i < nidx; i++) indices[i] = i; /* needs to be initialized */
14729566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(isnumb, &numbering));
14739566063dSJacob Faibussowitsch     PetscCall(PetscSortIntWithPermutation(nidx, numbering, indices));
14749566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(isnumb, &numbering));
147592bb6962SLisandro Dalcin     if (isbaij && bs > 1) { /* adjust for the block-aij case */
14769566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(nidx * bs, &newidx));
14772fa5cd67SKarl Rupp       for (i = 0; i < nidx; i++) {
14782fa5cd67SKarl Rupp         for (j = 0; j < bs; j++) newidx[i * bs + j] = indices[i] * bs + j;
14792fa5cd67SKarl Rupp       }
14809566063dSJacob Faibussowitsch       PetscCall(PetscFree(indices));
148192bb6962SLisandro Dalcin       nidx *= bs;
148292bb6962SLisandro Dalcin       indices = newidx;
148392bb6962SLisandro Dalcin     }
148492bb6962SLisandro Dalcin     /* Shift to get global indices */
148592bb6962SLisandro Dalcin     for (i = 0; i < nidx; i++) indices[i] += rstart;
148692bb6962SLisandro Dalcin 
148792bb6962SLisandro Dalcin     /* Build the index sets for each block */
148892bb6962SLisandro Dalcin     for (i = 0; i < n; i++) {
14899566063dSJacob Faibussowitsch       PetscCall(ISCreateGeneral(PETSC_COMM_SELF, count[i], &indices[start], PETSC_COPY_VALUES, &is[i]));
14909566063dSJacob Faibussowitsch       PetscCall(ISSort(is[i]));
149192bb6962SLisandro Dalcin       start += count[i];
149292bb6962SLisandro Dalcin     }
149392bb6962SLisandro Dalcin 
14949566063dSJacob Faibussowitsch     PetscCall(PetscFree(count));
14959566063dSJacob Faibussowitsch     PetscCall(PetscFree(indices));
14969566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&isnumb));
14979566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&ispart));
149892bb6962SLisandro Dalcin   }
14993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
150092bb6962SLisandro Dalcin }
150192bb6962SLisandro Dalcin 
150292bb6962SLisandro Dalcin /*@C
150392bb6962SLisandro Dalcin   PCASMDestroySubdomains - Destroys the index sets created with
1504f1580f4eSBarry Smith   `PCASMCreateSubdomains()`. Should be called after setting subdomains with `PCASMSetLocalSubdomains()`.
150592bb6962SLisandro Dalcin 
150692bb6962SLisandro Dalcin   Collective
150792bb6962SLisandro Dalcin 
150892bb6962SLisandro Dalcin   Input Parameters:
150992bb6962SLisandro Dalcin + n        - the number of index sets
15102b691e39SMatthew Knepley . is       - the array of index sets
15112fe279fdSBarry Smith - is_local - the array of local index sets, can be `NULL`
151292bb6962SLisandro Dalcin 
151392bb6962SLisandro Dalcin   Level: advanced
151492bb6962SLisandro Dalcin 
1515ce78bad3SBarry Smith   Developer Note:
1516ce78bad3SBarry Smith   The `IS` arguments should be a *[]
1517ce78bad3SBarry Smith 
1518562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMCreateSubdomains()`, `PCASMSetLocalSubdomains()`
151992bb6962SLisandro Dalcin @*/
1520ce78bad3SBarry Smith PetscErrorCode PCASMDestroySubdomains(PetscInt n, IS *is[], IS *is_local[])
1521d71ae5a4SJacob Faibussowitsch {
152292bb6962SLisandro Dalcin   PetscInt i;
15235fd66863SKarl Rupp 
152492bb6962SLisandro Dalcin   PetscFunctionBegin;
15253ba16761SJacob Faibussowitsch   if (n <= 0) PetscFunctionReturn(PETSC_SUCCESS);
1526ce78bad3SBarry Smith   if (*is) {
1527ce78bad3SBarry Smith     PetscAssertPointer(*is, 2);
1528ce78bad3SBarry Smith     for (i = 0; i < n; i++) PetscCall(ISDestroy(&(*is)[i]));
1529ce78bad3SBarry Smith     PetscCall(PetscFree(*is));
1530a35b7b57SDmitry Karpeev   }
1531ce78bad3SBarry Smith   if (is_local && *is_local) {
1532ce78bad3SBarry Smith     PetscAssertPointer(*is_local, 3);
1533ce78bad3SBarry Smith     for (i = 0; i < n; i++) PetscCall(ISDestroy(&(*is_local)[i]));
1534ce78bad3SBarry Smith     PetscCall(PetscFree(*is_local));
15352b691e39SMatthew Knepley   }
15363ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
153792bb6962SLisandro Dalcin }
153892bb6962SLisandro Dalcin 
15396141accfSBarry Smith /*@C
15404b9ad928SBarry Smith   PCASMCreateSubdomains2D - Creates the index sets for the overlapping Schwarz
1541f1580f4eSBarry Smith   preconditioner, `PCASM`, for a two-dimensional problem on a regular grid.
15424b9ad928SBarry Smith 
15434b9ad928SBarry Smith   Not Collective
15444b9ad928SBarry Smith 
15454b9ad928SBarry Smith   Input Parameters:
15466b867d5aSJose E. Roman + m       - the number of mesh points in the x direction
15476b867d5aSJose E. Roman . n       - the number of mesh points in the y direction
15486b867d5aSJose E. Roman . M       - the number of subdomains in the x direction
15496b867d5aSJose E. Roman . N       - the number of subdomains in the y direction
15504b9ad928SBarry Smith . dof     - degrees of freedom per node
15514b9ad928SBarry Smith - overlap - overlap in mesh lines
15524b9ad928SBarry Smith 
15534b9ad928SBarry Smith   Output Parameters:
15544b9ad928SBarry Smith + Nsub     - the number of subdomains created
15553d03bb48SJed Brown . is       - array of index sets defining overlapping (if overlap > 0) subdomains
15563d03bb48SJed Brown - is_local - array of index sets defining non-overlapping subdomains
15574b9ad928SBarry Smith 
155820f4b53cSBarry Smith   Level: advanced
155920f4b53cSBarry Smith 
15604b9ad928SBarry Smith   Note:
1561f1580f4eSBarry Smith   Presently `PCAMSCreateSubdomains2d()` is valid only for sequential
15624b9ad928SBarry Smith   preconditioners.  More general related routines are
1563f1580f4eSBarry Smith   `PCASMSetTotalSubdomains()` and `PCASMSetLocalSubdomains()`.
15644b9ad928SBarry Smith 
1565562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMSetTotalSubdomains()`, `PCASMSetLocalSubdomains()`, `PCASMGetSubKSP()`,
1566db781477SPatrick Sanan           `PCASMSetOverlap()`
15674b9ad928SBarry Smith @*/
1568ce78bad3SBarry Smith PetscErrorCode PCASMCreateSubdomains2D(PetscInt m, PetscInt n, PetscInt M, PetscInt N, PetscInt dof, PetscInt overlap, PetscInt *Nsub, IS *is[], IS *is_local[])
1569d71ae5a4SJacob Faibussowitsch {
15703d03bb48SJed Brown   PetscInt i, j, height, width, ystart, xstart, yleft, yright, xleft, xright, loc_outer;
157113f74950SBarry Smith   PetscInt nidx, *idx, loc, ii, jj, count;
15724b9ad928SBarry Smith 
15734b9ad928SBarry Smith   PetscFunctionBegin;
15747827d75bSBarry Smith   PetscCheck(dof == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "dof must be 1");
15754b9ad928SBarry Smith 
15764b9ad928SBarry Smith   *Nsub = N * M;
15779566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(*Nsub, is));
15789566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(*Nsub, is_local));
15794b9ad928SBarry Smith   ystart    = 0;
15803d03bb48SJed Brown   loc_outer = 0;
15814b9ad928SBarry Smith   for (i = 0; i < N; i++) {
15824b9ad928SBarry Smith     height = n / N + ((n % N) > i); /* height of subdomain */
15837827d75bSBarry Smith     PetscCheck(height >= 2, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many N subdomains for mesh dimension n");
15849371c9d4SSatish Balay     yleft = ystart - overlap;
15859371c9d4SSatish Balay     if (yleft < 0) yleft = 0;
15869371c9d4SSatish Balay     yright = ystart + height + overlap;
15879371c9d4SSatish Balay     if (yright > n) yright = n;
15884b9ad928SBarry Smith     xstart = 0;
15894b9ad928SBarry Smith     for (j = 0; j < M; j++) {
15904b9ad928SBarry Smith       width = m / M + ((m % M) > j); /* width of subdomain */
15917827d75bSBarry Smith       PetscCheck(width >= 2, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many M subdomains for mesh dimension m");
15929371c9d4SSatish Balay       xleft = xstart - overlap;
15939371c9d4SSatish Balay       if (xleft < 0) xleft = 0;
15949371c9d4SSatish Balay       xright = xstart + width + overlap;
15959371c9d4SSatish Balay       if (xright > m) xright = m;
15964b9ad928SBarry Smith       nidx = (xright - xleft) * (yright - yleft);
15979566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(nidx, &idx));
15984b9ad928SBarry Smith       loc = 0;
15994b9ad928SBarry Smith       for (ii = yleft; ii < yright; ii++) {
16004b9ad928SBarry Smith         count = m * ii + xleft;
16012fa5cd67SKarl Rupp         for (jj = xleft; jj < xright; jj++) idx[loc++] = count++;
16024b9ad928SBarry Smith       }
16039566063dSJacob Faibussowitsch       PetscCall(ISCreateGeneral(PETSC_COMM_SELF, nidx, idx, PETSC_COPY_VALUES, (*is) + loc_outer));
16043d03bb48SJed Brown       if (overlap == 0) {
16059566063dSJacob Faibussowitsch         PetscCall(PetscObjectReference((PetscObject)(*is)[loc_outer]));
16062fa5cd67SKarl Rupp 
16073d03bb48SJed Brown         (*is_local)[loc_outer] = (*is)[loc_outer];
16083d03bb48SJed Brown       } else {
16093d03bb48SJed Brown         for (loc = 0, ii = ystart; ii < ystart + height; ii++) {
1610ad540459SPierre Jolivet           for (jj = xstart; jj < xstart + width; jj++) idx[loc++] = m * ii + jj;
16113d03bb48SJed Brown         }
16129566063dSJacob Faibussowitsch         PetscCall(ISCreateGeneral(PETSC_COMM_SELF, loc, idx, PETSC_COPY_VALUES, *is_local + loc_outer));
16133d03bb48SJed Brown       }
16149566063dSJacob Faibussowitsch       PetscCall(PetscFree(idx));
16154b9ad928SBarry Smith       xstart += width;
16163d03bb48SJed Brown       loc_outer++;
16174b9ad928SBarry Smith     }
16184b9ad928SBarry Smith     ystart += height;
16194b9ad928SBarry Smith   }
16209566063dSJacob Faibussowitsch   for (i = 0; i < *Nsub; i++) PetscCall(ISSort((*is)[i]));
16213ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
16224b9ad928SBarry Smith }
16234b9ad928SBarry Smith 
16244b9ad928SBarry Smith /*@C
16254b9ad928SBarry Smith   PCASMGetLocalSubdomains - Gets the local subdomains (for this processor
1626f1580f4eSBarry Smith   only) for the additive Schwarz preconditioner, `PCASM`.
16274b9ad928SBarry Smith 
1628ad4df100SBarry Smith   Not Collective
16294b9ad928SBarry Smith 
16304b9ad928SBarry Smith   Input Parameter:
16314b9ad928SBarry Smith . pc - the preconditioner context
16324b9ad928SBarry Smith 
16334b9ad928SBarry Smith   Output Parameters:
1634f17a6784SPierre Jolivet + n        - if requested, the number of subdomains for this processor (default value = 1)
1635f17a6784SPierre Jolivet . is       - if requested, the index sets that define the subdomains for this processor
163620f4b53cSBarry Smith - is_local - if requested, the index sets that define the local part of the subdomains for this processor (can be `NULL`)
163720f4b53cSBarry Smith 
163820f4b53cSBarry Smith   Level: advanced
16394b9ad928SBarry Smith 
1640f1580f4eSBarry Smith   Note:
1641f1580f4eSBarry Smith   The `IS` numbering is in the parallel, global numbering of the vector.
16424b9ad928SBarry Smith 
1643562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMSetTotalSubdomains()`, `PCASMSetOverlap()`, `PCASMGetSubKSP()`,
1644db781477SPatrick Sanan           `PCASMCreateSubdomains2D()`, `PCASMSetLocalSubdomains()`, `PCASMGetLocalSubmatrices()`
16454b9ad928SBarry Smith @*/
1646d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMGetLocalSubdomains(PC pc, PetscInt *n, IS *is[], IS *is_local[])
1647d71ae5a4SJacob Faibussowitsch {
16482a808120SBarry Smith   PC_ASM   *osm = (PC_ASM *)pc->data;
1649ace3abfcSBarry Smith   PetscBool match;
16504b9ad928SBarry Smith 
16514b9ad928SBarry Smith   PetscFunctionBegin;
16520700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
16534f572ea9SToby Isaac   if (n) PetscAssertPointer(n, 2);
16544f572ea9SToby Isaac   if (is) PetscAssertPointer(is, 3);
16554f572ea9SToby Isaac   if (is_local) PetscAssertPointer(is_local, 4);
16569566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCASM, &match));
165728b400f6SJacob Faibussowitsch   PetscCheck(match, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "PC is not a PCASM");
16584b9ad928SBarry Smith   if (n) *n = osm->n_local_true;
16594b9ad928SBarry Smith   if (is) *is = osm->is;
16602b691e39SMatthew Knepley   if (is_local) *is_local = osm->is_local;
16613ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
16624b9ad928SBarry Smith }
16634b9ad928SBarry Smith 
16644b9ad928SBarry Smith /*@C
16654b9ad928SBarry Smith   PCASMGetLocalSubmatrices - Gets the local submatrices (for this processor
1666f1580f4eSBarry Smith   only) for the additive Schwarz preconditioner, `PCASM`.
16674b9ad928SBarry Smith 
1668ad4df100SBarry Smith   Not Collective
16694b9ad928SBarry Smith 
16704b9ad928SBarry Smith   Input Parameter:
16714b9ad928SBarry Smith . pc - the preconditioner context
16724b9ad928SBarry Smith 
16734b9ad928SBarry Smith   Output Parameters:
1674f17a6784SPierre Jolivet + n   - if requested, the number of matrices for this processor (default value = 1)
1675f17a6784SPierre Jolivet - mat - if requested, the matrices
16764b9ad928SBarry Smith 
16774b9ad928SBarry Smith   Level: advanced
16784b9ad928SBarry Smith 
167995452b02SPatrick Sanan   Notes:
1680f1580f4eSBarry Smith   Call after `PCSetUp()` (or `KSPSetUp()`) but before `PCApply()` and before `PCSetUpOnBlocks()`)
1681cf739d55SBarry Smith 
1682f1580f4eSBarry Smith   Usually one would use `PCSetModifySubMatrices()` to change the submatrices in building the preconditioner.
1683cf739d55SBarry Smith 
1684562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMSetTotalSubdomains()`, `PCASMSetOverlap()`, `PCASMGetSubKSP()`,
1685db781477SPatrick Sanan           `PCASMCreateSubdomains2D()`, `PCASMSetLocalSubdomains()`, `PCASMGetLocalSubdomains()`, `PCSetModifySubMatrices()`
16864b9ad928SBarry Smith @*/
1687d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMGetLocalSubmatrices(PC pc, PetscInt *n, Mat *mat[])
1688d71ae5a4SJacob Faibussowitsch {
16894b9ad928SBarry Smith   PC_ASM   *osm;
1690ace3abfcSBarry Smith   PetscBool match;
16914b9ad928SBarry Smith 
16924b9ad928SBarry Smith   PetscFunctionBegin;
16930700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
16944f572ea9SToby Isaac   if (n) PetscAssertPointer(n, 2);
16954f572ea9SToby Isaac   if (mat) PetscAssertPointer(mat, 3);
169628b400f6SJacob Faibussowitsch   PetscCheck(pc->setupcalled, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must call after KSPSetUp() or PCSetUp().");
16979566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCASM, &match));
169892bb6962SLisandro Dalcin   if (!match) {
169992bb6962SLisandro Dalcin     if (n) *n = 0;
17000298fd71SBarry Smith     if (mat) *mat = NULL;
170192bb6962SLisandro Dalcin   } else {
17024b9ad928SBarry Smith     osm = (PC_ASM *)pc->data;
17034b9ad928SBarry Smith     if (n) *n = osm->n_local_true;
17044b9ad928SBarry Smith     if (mat) *mat = osm->pmat;
170592bb6962SLisandro Dalcin   }
17063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
17074b9ad928SBarry Smith }
1708d709ea83SDmitry Karpeev 
1709d709ea83SDmitry Karpeev /*@
1710f1580f4eSBarry Smith   PCASMSetDMSubdomains - Indicates whether to use `DMCreateDomainDecomposition()` to define the subdomains, whenever possible.
1711f1ee410cSBarry Smith 
1712d709ea83SDmitry Karpeev   Logically Collective
1713d709ea83SDmitry Karpeev 
1714d8d19677SJose E. Roman   Input Parameters:
1715d709ea83SDmitry Karpeev + pc  - the preconditioner
1716f1580f4eSBarry Smith - flg - boolean indicating whether to use subdomains defined by the `DM`
1717d709ea83SDmitry Karpeev 
1718d709ea83SDmitry Karpeev   Options Database Key:
171973ff1848SBarry Smith . -pc_asm_dm_subdomains <bool> - use subdomains defined by the `DM` with `DMCreateDomainDecomposition()`
1720d709ea83SDmitry Karpeev 
1721d709ea83SDmitry Karpeev   Level: intermediate
1722d709ea83SDmitry Karpeev 
1723f1580f4eSBarry Smith   Note:
1724f1580f4eSBarry Smith   `PCASMSetTotalSubdomains()` and `PCASMSetOverlap()` take precedence over `PCASMSetDMSubdomains()`,
1725d709ea83SDmitry Karpeev   so setting either of the first two effectively turns the latter off.
1726d709ea83SDmitry Karpeev 
172773ff1848SBarry Smith   Developer Note:
172873ff1848SBarry Smith   This should be `PCASMSetUseDMSubdomains()`, similarly for the options database key
172973ff1848SBarry Smith 
1730562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMGetDMSubdomains()`, `PCASMSetTotalSubdomains()`, `PCASMSetOverlap()`
1731db781477SPatrick Sanan           `PCASMCreateSubdomains2D()`, `PCASMSetLocalSubdomains()`, `PCASMGetLocalSubdomains()`
1732d709ea83SDmitry Karpeev @*/
1733d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMSetDMSubdomains(PC pc, PetscBool flg)
1734d71ae5a4SJacob Faibussowitsch {
1735d709ea83SDmitry Karpeev   PC_ASM   *osm = (PC_ASM *)pc->data;
1736d709ea83SDmitry Karpeev   PetscBool match;
1737d709ea83SDmitry Karpeev 
1738d709ea83SDmitry Karpeev   PetscFunctionBegin;
1739d709ea83SDmitry Karpeev   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1740d709ea83SDmitry Karpeev   PetscValidLogicalCollectiveBool(pc, flg, 2);
174128b400f6SJacob Faibussowitsch   PetscCheck(!pc->setupcalled, ((PetscObject)pc)->comm, PETSC_ERR_ARG_WRONGSTATE, "Not for a setup PC.");
17429566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCASM, &match));
1743ad540459SPierre Jolivet   if (match) osm->dm_subdomains = flg;
17443ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1745d709ea83SDmitry Karpeev }
1746d709ea83SDmitry Karpeev 
1747d709ea83SDmitry Karpeev /*@
1748f1580f4eSBarry Smith   PCASMGetDMSubdomains - Returns flag indicating whether to use `DMCreateDomainDecomposition()` to define the subdomains, whenever possible.
1749f1580f4eSBarry Smith 
1750d709ea83SDmitry Karpeev   Not Collective
1751d709ea83SDmitry Karpeev 
1752d709ea83SDmitry Karpeev   Input Parameter:
1753d709ea83SDmitry Karpeev . pc - the preconditioner
1754d709ea83SDmitry Karpeev 
1755d709ea83SDmitry Karpeev   Output Parameter:
175620f4b53cSBarry Smith . flg - boolean indicating whether to use subdomains defined by the `DM`
1757d709ea83SDmitry Karpeev 
1758d709ea83SDmitry Karpeev   Level: intermediate
1759d709ea83SDmitry Karpeev 
176073ff1848SBarry Smith   Developer Note:
176173ff1848SBarry Smith   This should be `PCASMSetUseDMSubdomains()`
176273ff1848SBarry Smith 
1763562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMSetDMSubdomains()`, `PCASMSetTotalSubdomains()`, `PCASMSetOverlap()`
1764db781477SPatrick Sanan           `PCASMCreateSubdomains2D()`, `PCASMSetLocalSubdomains()`, `PCASMGetLocalSubdomains()`
1765d709ea83SDmitry Karpeev @*/
1766d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMGetDMSubdomains(PC pc, PetscBool *flg)
1767d71ae5a4SJacob Faibussowitsch {
1768d709ea83SDmitry Karpeev   PC_ASM   *osm = (PC_ASM *)pc->data;
1769d709ea83SDmitry Karpeev   PetscBool match;
1770d709ea83SDmitry Karpeev 
1771d709ea83SDmitry Karpeev   PetscFunctionBegin;
1772d709ea83SDmitry Karpeev   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
17734f572ea9SToby Isaac   PetscAssertPointer(flg, 2);
17749566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCASM, &match));
177556ea09ceSStefano Zampini   if (match) *flg = osm->dm_subdomains;
177656ea09ceSStefano Zampini   else *flg = PETSC_FALSE;
17773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1778d709ea83SDmitry Karpeev }
177980ec0b7dSPatrick Sanan 
17805d83a8b1SBarry Smith /*@
1781f1580f4eSBarry Smith   PCASMGetSubMatType - Gets the matrix type used for `PCASM` subsolves, as a string.
178280ec0b7dSPatrick Sanan 
178380ec0b7dSPatrick Sanan   Not Collective
178480ec0b7dSPatrick Sanan 
178580ec0b7dSPatrick Sanan   Input Parameter:
1786f1580f4eSBarry Smith . pc - the `PC`
178780ec0b7dSPatrick Sanan 
178880ec0b7dSPatrick Sanan   Output Parameter:
1789feefa0e1SJacob Faibussowitsch . sub_mat_type - name of matrix type
179080ec0b7dSPatrick Sanan 
179180ec0b7dSPatrick Sanan   Level: advanced
179280ec0b7dSPatrick Sanan 
1793562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMSetSubMatType()`, `PCSetType()`, `VecSetType()`, `MatType`, `Mat`
179480ec0b7dSPatrick Sanan @*/
1795d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMGetSubMatType(PC pc, MatType *sub_mat_type)
1796d71ae5a4SJacob Faibussowitsch {
179756ea09ceSStefano Zampini   PetscFunctionBegin;
179856ea09ceSStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1799cac4c232SBarry Smith   PetscTryMethod(pc, "PCASMGetSubMatType_C", (PC, MatType *), (pc, sub_mat_type));
18003ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
180180ec0b7dSPatrick Sanan }
180280ec0b7dSPatrick Sanan 
18035d83a8b1SBarry Smith /*@
1804f1580f4eSBarry Smith   PCASMSetSubMatType - Set the type of matrix used for `PCASM` subsolves
180580ec0b7dSPatrick Sanan 
1806c3339decSBarry Smith   Collective
180780ec0b7dSPatrick Sanan 
180880ec0b7dSPatrick Sanan   Input Parameters:
1809f1580f4eSBarry Smith + pc           - the `PC` object
1810f1580f4eSBarry Smith - sub_mat_type - the `MatType`
181180ec0b7dSPatrick Sanan 
181280ec0b7dSPatrick Sanan   Options Database Key:
1813f1580f4eSBarry Smith . -pc_asm_sub_mat_type  <sub_mat_type> - Sets the matrix type used for subsolves, for example, seqaijviennacl.
1814f1580f4eSBarry Smith    If you specify a base name like aijviennacl, the corresponding sequential type is assumed.
181580ec0b7dSPatrick Sanan 
1816f1580f4eSBarry Smith   Note:
18172fe279fdSBarry Smith   See `MatType` for available types
181880ec0b7dSPatrick Sanan 
181980ec0b7dSPatrick Sanan   Level: advanced
182080ec0b7dSPatrick Sanan 
1821562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMGetSubMatType()`, `PCSetType()`, `VecSetType()`, `MatType`, `Mat`
182280ec0b7dSPatrick Sanan @*/
1823d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMSetSubMatType(PC pc, MatType sub_mat_type)
1824d71ae5a4SJacob Faibussowitsch {
182556ea09ceSStefano Zampini   PetscFunctionBegin;
182656ea09ceSStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1827cac4c232SBarry Smith   PetscTryMethod(pc, "PCASMSetSubMatType_C", (PC, MatType), (pc, sub_mat_type));
18283ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
182980ec0b7dSPatrick Sanan }
1830