xref: /petsc/src/ksp/pc/impls/asm/asm.c (revision 3220ff8572602716d60bdda8b51773ebceb3c8ea)
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 
PCView_ASM(PC pc,PetscViewer viewer)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;
219f196a02SMartin Diehl   PetscBool         isascii, isstring;
224b9ad928SBarry Smith   PetscViewer       sviewer;
23ed155784SPierre Jolivet   PetscViewerFormat format;
24ed155784SPierre Jolivet   const char       *prefix;
254b9ad928SBarry Smith 
264b9ad928SBarry Smith   PetscFunctionBegin;
279f196a02SMartin Diehl   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
289566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSTRING, &isstring));
299f196a02SMartin Diehl   if (isascii) {
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 
PCASMPrintSubdomains(PC pc)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 
PCSetUp_ASM(PC pc)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]));
262*bf0c7fc2SBarry Smith           PetscCall(KSPSetDMActive(ksp, KSP_DMACTIVE_ALL, 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 
PCSetUpOnBlocks_ASM(PC pc)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 
PCApply_ASM(PC pc,Vec x,Vec y)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 
PCMatApply_ASM_Private(PC pc,Mat X,Mat Y,PetscBool transpose)4934dbf25a8SPierre 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   */
5074dbf25a8SPierre 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   }
5124dbf25a8SPierre 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 */
5374dbf25a8SPierre 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));
5414dbf25a8SPierre Jolivet   } else {
5424dbf25a8SPierre Jolivet     PetscCall(PetscLogEventBegin(PC_ApplyTransposeOnBlocks, osm->ksp[0], Z, W, 0));
5434dbf25a8SPierre Jolivet     PetscCall(KSPMatSolveTranspose(osm->ksp[0], Z, W));
5444dbf25a8SPierre Jolivet     PetscCall(PetscLogEventEnd(PC_ApplyTransposeOnBlocks, osm->ksp[0], Z, W, 0));
5454dbf25a8SPierre Jolivet   }
5464dbf25a8SPierre 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));
5524dbf25a8SPierre 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 
PCMatApply_ASM(PC pc,Mat X,Mat Y)5714dbf25a8SPierre Jolivet static PetscErrorCode PCMatApply_ASM(PC pc, Mat X, Mat Y)
5724dbf25a8SPierre Jolivet {
5734dbf25a8SPierre Jolivet   PetscFunctionBegin;
5744dbf25a8SPierre Jolivet   PetscCall(PCMatApply_ASM_Private(pc, X, Y, PETSC_FALSE));
5754dbf25a8SPierre Jolivet   PetscFunctionReturn(PETSC_SUCCESS);
5764dbf25a8SPierre Jolivet }
5774dbf25a8SPierre Jolivet 
PCMatApplyTranspose_ASM(PC pc,Mat X,Mat Y)5784dbf25a8SPierre Jolivet static PetscErrorCode PCMatApplyTranspose_ASM(PC pc, Mat X, Mat Y)
5794dbf25a8SPierre Jolivet {
5804dbf25a8SPierre Jolivet   PetscFunctionBegin;
5814dbf25a8SPierre Jolivet   PetscCall(PCMatApply_ASM_Private(pc, X, Y, PETSC_TRUE));
5824dbf25a8SPierre Jolivet   PetscFunctionReturn(PETSC_SUCCESS);
5834dbf25a8SPierre Jolivet }
5844dbf25a8SPierre Jolivet 
PCApplyTranspose_ASM(PC pc,Vec x,Vec y)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;
5926040bc1fSPierre Jolivet   PetscCheck(osm->n_local_true <= 1 || osm->loctype == PC_COMPOSITE_ADDITIVE, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Not yet implemented");
5934b9ad928SBarry Smith   /*
5944b9ad928SBarry Smith      Support for limiting the restriction or interpolation to only local
5954b9ad928SBarry Smith      subdomain values (leaving the other values 0).
5964b9ad928SBarry Smith 
5974b9ad928SBarry Smith      Note: these are reversed from the PCApply_ASM() because we are applying the
5984b9ad928SBarry Smith      transpose of the three terms
5994b9ad928SBarry Smith   */
600d8b3b5e3Seaulisa 
6014b9ad928SBarry Smith   if (!(osm->type & PC_ASM_INTERPOLATE)) {
6024b9ad928SBarry Smith     forward = SCATTER_FORWARD_LOCAL;
6034b9ad928SBarry Smith     /* have to zero the work RHS since scatter may leave some slots empty */
6049566063dSJacob Faibussowitsch     PetscCall(VecSet(osm->lx, 0.0));
6054b9ad928SBarry Smith   }
6062fa5cd67SKarl Rupp   if (!(osm->type & PC_ASM_RESTRICT)) reverse = SCATTER_REVERSE_LOCAL;
6074b9ad928SBarry Smith 
608b0de9f23SBarry Smith   /* zero the global and the local solutions */
6099566063dSJacob Faibussowitsch   PetscCall(VecSet(y, 0.0));
6109566063dSJacob Faibussowitsch   PetscCall(VecSet(osm->ly, 0.0));
611d8b3b5e3Seaulisa 
612b0de9f23SBarry Smith   /* Copy the global RHS to local RHS including the ghost nodes */
6139566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(osm->restriction, x, osm->lx, INSERT_VALUES, forward));
6149566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(osm->restriction, x, osm->lx, INSERT_VALUES, forward));
615d8b3b5e3Seaulisa 
616b0de9f23SBarry Smith   /* Restrict local RHS to the overlapping 0-block RHS */
6179566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward));
6189566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward));
619d8b3b5e3Seaulisa 
6204b9ad928SBarry Smith   /* do the local solves */
621d8b3b5e3Seaulisa   for (i = 0; i < n_local_true; ++i) {
622b0de9f23SBarry Smith     /* solve the overlapping i-block */
6236040bc1fSPierre Jolivet     PetscCall(PetscLogEventBegin(PC_ApplyTransposeOnBlocks, osm->ksp[i], osm->x[i], osm->y[i], 0));
6249566063dSJacob Faibussowitsch     PetscCall(KSPSolveTranspose(osm->ksp[i], osm->x[i], osm->y[i]));
6259566063dSJacob Faibussowitsch     PetscCall(KSPCheckSolve(osm->ksp[i], pc, osm->y[i]));
6266040bc1fSPierre Jolivet     PetscCall(PetscLogEventEnd(PC_ApplyTransposeOnBlocks, osm->ksp[i], osm->x[i], osm->y[i], 0));
627d8b3b5e3Seaulisa 
628a681a0f1SPierre Jolivet     if (osm->lprolongation && osm->type != PC_ASM_RESTRICT) { /* interpolate the non-overlapping i-block solution to the local solution */
6299566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward));
6309566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward));
631910cf402Sprj-     } else { /* interpolate the overlapping i-block solution to the local solution */
6329566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse));
6339566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse));
6344b9ad928SBarry Smith     }
635d8b3b5e3Seaulisa 
636d8b3b5e3Seaulisa     if (i < n_local_true - 1) {
637b0de9f23SBarry Smith       /* Restrict local RHS to the overlapping (i+1)-block RHS */
6389566063dSJacob Faibussowitsch       PetscCall(VecScatterBegin(osm->lrestriction[i + 1], osm->lx, osm->x[i + 1], INSERT_VALUES, forward));
6399566063dSJacob Faibussowitsch       PetscCall(VecScatterEnd(osm->lrestriction[i + 1], osm->lx, osm->x[i + 1], INSERT_VALUES, forward));
6404b9ad928SBarry Smith     }
6414b9ad928SBarry Smith   }
642b0de9f23SBarry Smith   /* Add the local solution to the global solution including the ghost nodes */
6439566063dSJacob Faibussowitsch   PetscCall(VecScatterBegin(osm->restriction, osm->ly, y, ADD_VALUES, reverse));
6449566063dSJacob Faibussowitsch   PetscCall(VecScatterEnd(osm->restriction, osm->ly, y, ADD_VALUES, reverse));
6453ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
6464b9ad928SBarry Smith }
6474b9ad928SBarry Smith 
PCReset_ASM(PC pc)648d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCReset_ASM(PC pc)
649d71ae5a4SJacob Faibussowitsch {
6504b9ad928SBarry Smith   PC_ASM  *osm = (PC_ASM *)pc->data;
65113f74950SBarry Smith   PetscInt i;
6524b9ad928SBarry Smith 
6534b9ad928SBarry Smith   PetscFunctionBegin;
65492bb6962SLisandro Dalcin   if (osm->ksp) {
65548a46eb9SPierre Jolivet     for (i = 0; i < osm->n_local_true; i++) PetscCall(KSPReset(osm->ksp[i]));
65692bb6962SLisandro Dalcin   }
657e09e08ccSBarry Smith   if (osm->pmat) {
65848a46eb9SPierre Jolivet     if (osm->n_local_true > 0) PetscCall(MatDestroySubMatrices(osm->n_local_true, &osm->pmat));
65992bb6962SLisandro Dalcin   }
6601dd8081eSeaulisa   if (osm->lrestriction) {
6619566063dSJacob Faibussowitsch     PetscCall(VecScatterDestroy(&osm->restriction));
6621dd8081eSeaulisa     for (i = 0; i < osm->n_local_true; i++) {
6639566063dSJacob Faibussowitsch       PetscCall(VecScatterDestroy(&osm->lrestriction[i]));
6649566063dSJacob Faibussowitsch       if (osm->lprolongation) PetscCall(VecScatterDestroy(&osm->lprolongation[i]));
6659566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&osm->x[i]));
6669566063dSJacob Faibussowitsch       PetscCall(VecDestroy(&osm->y[i]));
6674b9ad928SBarry Smith     }
6689566063dSJacob Faibussowitsch     PetscCall(PetscFree(osm->lrestriction));
6699566063dSJacob Faibussowitsch     if (osm->lprolongation) PetscCall(PetscFree(osm->lprolongation));
6709566063dSJacob Faibussowitsch     PetscCall(PetscFree(osm->x));
6719566063dSJacob Faibussowitsch     PetscCall(PetscFree(osm->y));
67292bb6962SLisandro Dalcin   }
673ce78bad3SBarry Smith   PetscCall(PCASMDestroySubdomains(osm->n_local_true, &osm->is, &osm->is_local));
6749566063dSJacob Faibussowitsch   PetscCall(ISDestroy(&osm->lis));
6759566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&osm->lx));
6769566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&osm->ly));
67748a46eb9SPierre Jolivet   if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE) PetscCall(MatDestroyMatrices(osm->n_local_true, &osm->lmats));
6782fa5cd67SKarl Rupp 
6799566063dSJacob Faibussowitsch   PetscCall(PetscFree(osm->sub_mat_type));
68080ec0b7dSPatrick Sanan 
6810a545947SLisandro Dalcin   osm->is       = NULL;
6820a545947SLisandro Dalcin   osm->is_local = NULL;
6833ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
684e91c6855SBarry Smith }
685e91c6855SBarry Smith 
PCDestroy_ASM(PC pc)686d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_ASM(PC pc)
687d71ae5a4SJacob Faibussowitsch {
688e91c6855SBarry Smith   PC_ASM  *osm = (PC_ASM *)pc->data;
689e91c6855SBarry Smith   PetscInt i;
690e91c6855SBarry Smith 
691e91c6855SBarry Smith   PetscFunctionBegin;
6929566063dSJacob Faibussowitsch   PetscCall(PCReset_ASM(pc));
693e91c6855SBarry Smith   if (osm->ksp) {
69448a46eb9SPierre Jolivet     for (i = 0; i < osm->n_local_true; i++) PetscCall(KSPDestroy(&osm->ksp[i]));
6959566063dSJacob Faibussowitsch     PetscCall(PetscFree(osm->ksp));
696e91c6855SBarry Smith   }
6979566063dSJacob Faibussowitsch   PetscCall(PetscFree(pc->data));
69896322394SPierre Jolivet 
6999566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetLocalSubdomains_C", NULL));
7009566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetTotalSubdomains_C", NULL));
7019566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetOverlap_C", NULL));
7029566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetType_C", NULL));
7039566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMGetType_C", NULL));
7049566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetLocalType_C", NULL));
7059566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMGetLocalType_C", NULL));
7069566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetSortIndices_C", NULL));
7079566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMGetSubKSP_C", NULL));
7089566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMGetSubMatType_C", NULL));
7099566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetSubMatType_C", NULL));
7103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7114b9ad928SBarry Smith }
7124b9ad928SBarry Smith 
PCSetFromOptions_ASM(PC pc,PetscOptionItems PetscOptionsObject)713ce78bad3SBarry Smith static PetscErrorCode PCSetFromOptions_ASM(PC pc, PetscOptionItems PetscOptionsObject)
714d71ae5a4SJacob Faibussowitsch {
7154b9ad928SBarry Smith   PC_ASM         *osm = (PC_ASM *)pc->data;
7169dcbbd2bSBarry Smith   PetscInt        blocks, ovl;
71757501b6eSBarry Smith   PetscBool       flg;
71892bb6962SLisandro Dalcin   PCASMType       asmtype;
71912cd4985SMatthew G. Knepley   PCCompositeType loctype;
72080ec0b7dSPatrick Sanan   char            sub_mat_type[256];
7214b9ad928SBarry Smith 
7224b9ad928SBarry Smith   PetscFunctionBegin;
723d0609cedSBarry Smith   PetscOptionsHeadBegin(PetscOptionsObject, "Additive Schwarz options");
7249566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-pc_asm_dm_subdomains", "Use DMCreateDomainDecomposition() to define subdomains", "PCASMSetDMSubdomains", osm->dm_subdomains, &osm->dm_subdomains, &flg));
7259566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_asm_blocks", "Number of subdomains", "PCASMSetTotalSubdomains", osm->n, &blocks, &flg));
72665db9045SDmitry Karpeev   if (flg) {
7279566063dSJacob Faibussowitsch     PetscCall(PCASMSetTotalSubdomains(pc, blocks, NULL, NULL));
728d709ea83SDmitry Karpeev     osm->dm_subdomains = PETSC_FALSE;
72965db9045SDmitry Karpeev   }
7309566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_asm_local_blocks", "Number of local subdomains", "PCASMSetLocalSubdomains", osm->n_local_true, &blocks, &flg));
731342c94f9SMatthew G. Knepley   if (flg) {
7329566063dSJacob Faibussowitsch     PetscCall(PCASMSetLocalSubdomains(pc, blocks, NULL, NULL));
733342c94f9SMatthew G. Knepley     osm->dm_subdomains = PETSC_FALSE;
734342c94f9SMatthew G. Knepley   }
7359566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-pc_asm_overlap", "Number of grid points overlap", "PCASMSetOverlap", osm->overlap, &ovl, &flg));
73665db9045SDmitry Karpeev   if (flg) {
7379566063dSJacob Faibussowitsch     PetscCall(PCASMSetOverlap(pc, ovl));
738d709ea83SDmitry Karpeev     osm->dm_subdomains = PETSC_FALSE;
73965db9045SDmitry Karpeev   }
74090d69ab7SBarry Smith   flg = PETSC_FALSE;
7419566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-pc_asm_type", "Type of restriction/extension", "PCASMSetType", PCASMTypes, (PetscEnum)osm->type, (PetscEnum *)&asmtype, &flg));
7429566063dSJacob Faibussowitsch   if (flg) PetscCall(PCASMSetType(pc, asmtype));
74312cd4985SMatthew G. Knepley   flg = PETSC_FALSE;
7449566063dSJacob Faibussowitsch   PetscCall(PetscOptionsEnum("-pc_asm_local_type", "Type of local solver composition", "PCASMSetLocalType", PCCompositeTypes, (PetscEnum)osm->loctype, (PetscEnum *)&loctype, &flg));
7459566063dSJacob Faibussowitsch   if (flg) PetscCall(PCASMSetLocalType(pc, loctype));
7469566063dSJacob Faibussowitsch   PetscCall(PetscOptionsFList("-pc_asm_sub_mat_type", "Subsolve Matrix Type", "PCASMSetSubMatType", MatList, NULL, sub_mat_type, 256, &flg));
7471baa6e33SBarry Smith   if (flg) PetscCall(PCASMSetSubMatType(pc, sub_mat_type));
748d0609cedSBarry Smith   PetscOptionsHeadEnd();
7493ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
7504b9ad928SBarry Smith }
7514b9ad928SBarry Smith 
PCASMSetLocalSubdomains_ASM(PC pc,PetscInt n,IS is[],IS is_local[])752d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCASMSetLocalSubdomains_ASM(PC pc, PetscInt n, IS is[], IS is_local[])
753d71ae5a4SJacob Faibussowitsch {
7544b9ad928SBarry Smith   PC_ASM  *osm = (PC_ASM *)pc->data;
75592bb6962SLisandro Dalcin   PetscInt i;
7564b9ad928SBarry Smith 
7574b9ad928SBarry Smith   PetscFunctionBegin;
75863a3b9bcSJacob Faibussowitsch   PetscCheck(n >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Each process must have 1 or more blocks, n = %" PetscInt_FMT, n);
7597827d75bSBarry Smith   PetscCheck(!pc->setupcalled || (n == osm->n_local_true && !is), PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "PCASMSetLocalSubdomains() should be called before calling PCSetUp().");
760e7e72b3dSBarry Smith 
7614b9ad928SBarry Smith   if (!pc->setupcalled) {
76292bb6962SLisandro Dalcin     if (is) {
7639566063dSJacob Faibussowitsch       for (i = 0; i < n; i++) PetscCall(PetscObjectReference((PetscObject)is[i]));
76492bb6962SLisandro Dalcin     }
765832fc9a5SMatthew Knepley     if (is_local) {
7669566063dSJacob Faibussowitsch       for (i = 0; i < n; i++) PetscCall(PetscObjectReference((PetscObject)is_local[i]));
767832fc9a5SMatthew Knepley     }
768ce78bad3SBarry Smith     PetscCall(PCASMDestroySubdomains(osm->n_local_true, &osm->is, &osm->is_local));
7692fa5cd67SKarl Rupp 
77007517c86SMark Adams     if (osm->ksp && osm->n_local_true != n) {
77107517c86SMark Adams       for (i = 0; i < osm->n_local_true; i++) PetscCall(KSPDestroy(&osm->ksp[i]));
77207517c86SMark Adams       PetscCall(PetscFree(osm->ksp));
77307517c86SMark Adams     }
77407517c86SMark Adams 
7754b9ad928SBarry Smith     osm->n_local_true = n;
7760a545947SLisandro Dalcin     osm->is           = NULL;
7770a545947SLisandro Dalcin     osm->is_local     = NULL;
77892bb6962SLisandro Dalcin     if (is) {
7799566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(n, &osm->is));
7802fa5cd67SKarl Rupp       for (i = 0; i < n; i++) osm->is[i] = is[i];
7813d03bb48SJed Brown       /* Flag indicating that the user has set overlapping subdomains so PCASM should not increase their size. */
7823d03bb48SJed Brown       osm->overlap = -1;
78392bb6962SLisandro Dalcin     }
7842b691e39SMatthew Knepley     if (is_local) {
7859566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(n, &osm->is_local));
7862fa5cd67SKarl Rupp       for (i = 0; i < n; i++) osm->is_local[i] = is_local[i];
787a35b7b57SDmitry Karpeev       if (!is) {
7889566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(osm->n_local_true, &osm->is));
789a35b7b57SDmitry Karpeev         for (i = 0; i < osm->n_local_true; i++) {
790a35b7b57SDmitry Karpeev           if (osm->overlap > 0) { /* With positive overlap, osm->is[i] will be modified */
7919566063dSJacob Faibussowitsch             PetscCall(ISDuplicate(osm->is_local[i], &osm->is[i]));
7929566063dSJacob Faibussowitsch             PetscCall(ISCopy(osm->is_local[i], osm->is[i]));
793a35b7b57SDmitry Karpeev           } else {
7949566063dSJacob Faibussowitsch             PetscCall(PetscObjectReference((PetscObject)osm->is_local[i]));
795a35b7b57SDmitry Karpeev             osm->is[i] = osm->is_local[i];
796a35b7b57SDmitry Karpeev           }
797a35b7b57SDmitry Karpeev         }
798a35b7b57SDmitry Karpeev       }
7992b691e39SMatthew Knepley     }
8004b9ad928SBarry Smith   }
8013ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8024b9ad928SBarry Smith }
8034b9ad928SBarry Smith 
PCASMSetTotalSubdomains_ASM(PC pc,PetscInt N,IS * is,IS * is_local)804d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCASMSetTotalSubdomains_ASM(PC pc, PetscInt N, IS *is, IS *is_local)
805d71ae5a4SJacob Faibussowitsch {
8064b9ad928SBarry Smith   PC_ASM     *osm = (PC_ASM *)pc->data;
80713f74950SBarry Smith   PetscMPIInt rank, size;
80878904715SLisandro Dalcin   PetscInt    n;
8094b9ad928SBarry Smith 
8104b9ad928SBarry Smith   PetscFunctionBegin;
81163a3b9bcSJacob Faibussowitsch   PetscCheck(N >= 1, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Number of total blocks must be > 0, N = %" PetscInt_FMT, N);
81200045ab3SPierre Jolivet   PetscCheck(!is && !is_local, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Use PCASMSetLocalSubdomains() to set specific index sets, they cannot be set globally yet.");
8134b9ad928SBarry Smith 
8144b9ad928SBarry Smith   /*
815880770e9SJed Brown      Split the subdomains equally among all processors
8164b9ad928SBarry Smith   */
8179566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank));
8189566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size));
8194b9ad928SBarry Smith   n = N / size + ((N % size) > rank);
820835f2295SStefano 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);
8217827d75bSBarry Smith   PetscCheck(!pc->setupcalled || n == osm->n_local_true, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "PCASMSetTotalSubdomains() should be called before PCSetUp().");
8224b9ad928SBarry Smith   if (!pc->setupcalled) {
823ce78bad3SBarry Smith     PetscCall(PCASMDestroySubdomains(osm->n_local_true, &osm->is, &osm->is_local));
8242fa5cd67SKarl Rupp 
8254b9ad928SBarry Smith     osm->n_local_true = n;
8260a545947SLisandro Dalcin     osm->is           = NULL;
8270a545947SLisandro Dalcin     osm->is_local     = NULL;
8284b9ad928SBarry Smith   }
8293ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8304b9ad928SBarry Smith }
8314b9ad928SBarry Smith 
PCASMSetOverlap_ASM(PC pc,PetscInt ovl)832d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCASMSetOverlap_ASM(PC pc, PetscInt ovl)
833d71ae5a4SJacob Faibussowitsch {
83492bb6962SLisandro Dalcin   PC_ASM *osm = (PC_ASM *)pc->data;
8354b9ad928SBarry Smith 
8364b9ad928SBarry Smith   PetscFunctionBegin;
8377827d75bSBarry Smith   PetscCheck(ovl >= 0, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Negative overlap value requested");
8387827d75bSBarry Smith   PetscCheck(!pc->setupcalled || ovl == osm->overlap, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "PCASMSetOverlap() should be called before PCSetUp().");
8392fa5cd67SKarl Rupp   if (!pc->setupcalled) osm->overlap = ovl;
8403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8414b9ad928SBarry Smith }
8424b9ad928SBarry Smith 
PCASMSetType_ASM(PC pc,PCASMType type)843d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCASMSetType_ASM(PC pc, PCASMType type)
844d71ae5a4SJacob Faibussowitsch {
84592bb6962SLisandro Dalcin   PC_ASM *osm = (PC_ASM *)pc->data;
8464b9ad928SBarry Smith 
8474b9ad928SBarry Smith   PetscFunctionBegin;
8484b9ad928SBarry Smith   osm->type     = type;
849bf108f30SBarry Smith   osm->type_set = PETSC_TRUE;
8503ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8514b9ad928SBarry Smith }
8524b9ad928SBarry Smith 
PCASMGetType_ASM(PC pc,PCASMType * type)853d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCASMGetType_ASM(PC pc, PCASMType *type)
854d71ae5a4SJacob Faibussowitsch {
855c60c7ad4SBarry Smith   PC_ASM *osm = (PC_ASM *)pc->data;
856c60c7ad4SBarry Smith 
857c60c7ad4SBarry Smith   PetscFunctionBegin;
858c60c7ad4SBarry Smith   *type = osm->type;
8593ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
860c60c7ad4SBarry Smith }
861c60c7ad4SBarry Smith 
PCASMSetLocalType_ASM(PC pc,PCCompositeType type)862d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCASMSetLocalType_ASM(PC pc, PCCompositeType type)
863d71ae5a4SJacob Faibussowitsch {
86412cd4985SMatthew G. Knepley   PC_ASM *osm = (PC_ASM *)pc->data;
86512cd4985SMatthew G. Knepley 
86612cd4985SMatthew G. Knepley   PetscFunctionBegin;
8677827d75bSBarry 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");
86812cd4985SMatthew G. Knepley   osm->loctype = type;
8693ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
87012cd4985SMatthew G. Knepley }
87112cd4985SMatthew G. Knepley 
PCASMGetLocalType_ASM(PC pc,PCCompositeType * type)872d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCASMGetLocalType_ASM(PC pc, PCCompositeType *type)
873d71ae5a4SJacob Faibussowitsch {
87412cd4985SMatthew G. Knepley   PC_ASM *osm = (PC_ASM *)pc->data;
87512cd4985SMatthew G. Knepley 
87612cd4985SMatthew G. Knepley   PetscFunctionBegin;
87712cd4985SMatthew G. Knepley   *type = osm->loctype;
8783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
87912cd4985SMatthew G. Knepley }
88012cd4985SMatthew G. Knepley 
PCASMSetSortIndices_ASM(PC pc,PetscBool doSort)881d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCASMSetSortIndices_ASM(PC pc, PetscBool doSort)
882d71ae5a4SJacob Faibussowitsch {
8836ed231c7SMatthew Knepley   PC_ASM *osm = (PC_ASM *)pc->data;
8846ed231c7SMatthew Knepley 
8856ed231c7SMatthew Knepley   PetscFunctionBegin;
8866ed231c7SMatthew Knepley   osm->sort_indices = doSort;
8873ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
8886ed231c7SMatthew Knepley }
8896ed231c7SMatthew Knepley 
PCASMGetSubKSP_ASM(PC pc,PetscInt * n_local,PetscInt * first_local,KSP ** ksp)890d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCASMGetSubKSP_ASM(PC pc, PetscInt *n_local, PetscInt *first_local, KSP **ksp)
891d71ae5a4SJacob Faibussowitsch {
89292bb6962SLisandro Dalcin   PC_ASM *osm = (PC_ASM *)pc->data;
8934b9ad928SBarry Smith 
8944b9ad928SBarry Smith   PetscFunctionBegin;
8957827d75bSBarry 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");
8964b9ad928SBarry Smith 
8972fa5cd67SKarl Rupp   if (n_local) *n_local = osm->n_local_true;
89892bb6962SLisandro Dalcin   if (first_local) {
8999566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Scan(&osm->n_local_true, first_local, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)pc)));
90092bb6962SLisandro Dalcin     *first_local -= osm->n_local_true;
90192bb6962SLisandro Dalcin   }
902ed155784SPierre Jolivet   if (ksp) *ksp = osm->ksp;
9033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9044b9ad928SBarry Smith }
9054b9ad928SBarry Smith 
PCASMGetSubMatType_ASM(PC pc,MatType * sub_mat_type)906d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCASMGetSubMatType_ASM(PC pc, MatType *sub_mat_type)
907d71ae5a4SJacob Faibussowitsch {
90880ec0b7dSPatrick Sanan   PC_ASM *osm = (PC_ASM *)pc->data;
90980ec0b7dSPatrick Sanan 
91080ec0b7dSPatrick Sanan   PetscFunctionBegin;
91180ec0b7dSPatrick Sanan   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
9124f572ea9SToby Isaac   PetscAssertPointer(sub_mat_type, 2);
91380ec0b7dSPatrick Sanan   *sub_mat_type = osm->sub_mat_type;
9143ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
91580ec0b7dSPatrick Sanan }
91680ec0b7dSPatrick Sanan 
PCASMSetSubMatType_ASM(PC pc,MatType sub_mat_type)917d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCASMSetSubMatType_ASM(PC pc, MatType sub_mat_type)
918d71ae5a4SJacob Faibussowitsch {
91980ec0b7dSPatrick Sanan   PC_ASM *osm = (PC_ASM *)pc->data;
92080ec0b7dSPatrick Sanan 
92180ec0b7dSPatrick Sanan   PetscFunctionBegin;
92280ec0b7dSPatrick Sanan   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
9239566063dSJacob Faibussowitsch   PetscCall(PetscFree(osm->sub_mat_type));
9249566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(sub_mat_type, (char **)&osm->sub_mat_type));
9253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
92680ec0b7dSPatrick Sanan }
92780ec0b7dSPatrick Sanan 
9285d83a8b1SBarry Smith /*@
929f1580f4eSBarry Smith   PCASMSetLocalSubdomains - Sets the local subdomains (for this processor only) for the additive Schwarz preconditioner `PCASM`.
9304b9ad928SBarry Smith 
931c3339decSBarry Smith   Collective
9324b9ad928SBarry Smith 
9334b9ad928SBarry Smith   Input Parameters:
9344b9ad928SBarry Smith + pc       - the preconditioner context
9354b9ad928SBarry Smith . n        - the number of subdomains for this processor (default value = 1)
936a3b724e8SBarry Smith . is       - the index set that defines the subdomains for this processor (or `NULL` for PETSc to determine subdomains)
937f13dfd9eSBarry Smith              the values of the `is` array are copied so you can free the array (not the `IS` in the array) after this call
938f13dfd9eSBarry 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`
939f13dfd9eSBarry Smith              (or `NULL` to not provide these). The values of the `is_local` array are copied so you can free the array
940f13dfd9eSBarry Smith              (not the `IS` in the array) after this call
9414b9ad928SBarry Smith 
942342c94f9SMatthew G. Knepley   Options Database Key:
94320f4b53cSBarry Smith . -pc_asm_local_blocks <blks> - Sets number of local blocks
94420f4b53cSBarry Smith 
94520f4b53cSBarry Smith   Level: advanced
946342c94f9SMatthew G. Knepley 
9474b9ad928SBarry Smith   Notes:
948a3b724e8SBarry Smith   The `IS` numbering is in the parallel, global numbering of the vector for both `is` and `is_local`
9494b9ad928SBarry Smith 
950f1580f4eSBarry Smith   By default the `PCASM` preconditioner uses 1 block per processor.
9514b9ad928SBarry Smith 
952f1580f4eSBarry Smith   Use `PCASMSetTotalSubdomains()` to set the subdomains for all processors.
9534b9ad928SBarry Smith 
954a3b724e8SBarry Smith   If `is_local` is provided and `PCASMType` is `PC_ASM_RESTRICT` then the solution only over the `is_local` region is interpolated
955a3b724e8SBarry Smith   back to form the global solution (this is the standard restricted additive Schwarz method, RASM)
956f1ee410cSBarry Smith 
957a3b724e8SBarry Smith   If `is_local` is provided and `PCASMType` is `PC_ASM_INTERPOLATE` or `PC_ASM_NONE` then an error is generated since there is
958f1ee410cSBarry Smith   no code to handle that case.
959f1ee410cSBarry Smith 
960562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMSetTotalSubdomains()`, `PCASMSetOverlap()`, `PCASMGetSubKSP()`,
961f1580f4eSBarry Smith           `PCASMCreateSubdomains2D()`, `PCASMGetLocalSubdomains()`, `PCASMType`, `PCASMSetType()`, `PCGASM`
9624b9ad928SBarry Smith @*/
PCASMSetLocalSubdomains(PC pc,PetscInt n,IS is[],IS is_local[])963d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMSetLocalSubdomains(PC pc, PetscInt n, IS is[], IS is_local[])
964d71ae5a4SJacob Faibussowitsch {
9654b9ad928SBarry Smith   PetscFunctionBegin;
9660700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
967cac4c232SBarry Smith   PetscTryMethod(pc, "PCASMSetLocalSubdomains_C", (PC, PetscInt, IS[], IS[]), (pc, n, is, is_local));
9683ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
9694b9ad928SBarry Smith }
9704b9ad928SBarry Smith 
9715d83a8b1SBarry Smith /*@
972feb221f8SDmitry Karpeev   PCASMSetTotalSubdomains - Sets the subdomains for all processors for the
973f1580f4eSBarry Smith   additive Schwarz preconditioner, `PCASM`.
9744b9ad928SBarry Smith 
975c3339decSBarry Smith   Collective, all MPI ranks must pass in the same array of `IS`
9764b9ad928SBarry Smith 
9774b9ad928SBarry Smith   Input Parameters:
9784b9ad928SBarry Smith + pc       - the preconditioner context
979feb221f8SDmitry Karpeev . N        - the number of subdomains for all processors
980a3b724e8SBarry Smith . is       - the index sets that define the subdomains for all processors (or `NULL` to ask PETSc to determine the subdomains)
9815d83a8b1SBarry Smith              the values of the `is` array are copied so you can free the array (not the `IS` in the array) after this call
982a3b724e8SBarry Smith - is_local - the index sets that define the local part of the subdomains for this processor (or `NULL` to not provide this information)
9835d83a8b1SBarry 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
9844b9ad928SBarry Smith 
9854b9ad928SBarry Smith   Options Database Key:
9864b9ad928SBarry Smith . -pc_asm_blocks <blks> - Sets total blocks
9874b9ad928SBarry Smith 
98820f4b53cSBarry Smith   Level: advanced
98920f4b53cSBarry Smith 
9904b9ad928SBarry Smith   Notes:
991a3b724e8SBarry Smith   Currently you cannot use this to set the actual subdomains with the argument `is` or `is_local`.
9924b9ad928SBarry Smith 
993f1580f4eSBarry Smith   By default the `PCASM` preconditioner uses 1 block per processor.
9944b9ad928SBarry Smith 
9954b9ad928SBarry Smith   These index sets cannot be destroyed until after completion of the
996f1580f4eSBarry Smith   linear solves for which the `PCASM` preconditioner is being used.
9974b9ad928SBarry Smith 
998f1580f4eSBarry Smith   Use `PCASMSetLocalSubdomains()` to set local subdomains.
9994b9ad928SBarry Smith 
1000f1580f4eSBarry Smith   The `IS` numbering is in the parallel, global numbering of the vector for both is and is_local
10011093a601SBarry Smith 
1002562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMSetLocalSubdomains()`, `PCASMSetOverlap()`, `PCASMGetSubKSP()`,
1003f1580f4eSBarry Smith           `PCASMCreateSubdomains2D()`, `PCGASM`
10044b9ad928SBarry Smith @*/
PCASMSetTotalSubdomains(PC pc,PetscInt N,IS is[],IS is_local[])1005d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMSetTotalSubdomains(PC pc, PetscInt N, IS is[], IS is_local[])
1006d71ae5a4SJacob Faibussowitsch {
10074b9ad928SBarry Smith   PetscFunctionBegin;
10080700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1009cac4c232SBarry Smith   PetscTryMethod(pc, "PCASMSetTotalSubdomains_C", (PC, PetscInt, IS[], IS[]), (pc, N, is, is_local));
10103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10114b9ad928SBarry Smith }
10124b9ad928SBarry Smith 
10134b9ad928SBarry Smith /*@
10144b9ad928SBarry Smith   PCASMSetOverlap - Sets the overlap between a pair of subdomains for the
1015f1580f4eSBarry Smith   additive Schwarz preconditioner, `PCASM`.
10164b9ad928SBarry Smith 
1017c3339decSBarry Smith   Logically Collective
10184b9ad928SBarry Smith 
10194b9ad928SBarry Smith   Input Parameters:
10204b9ad928SBarry Smith + pc  - the preconditioner context
10214b9ad928SBarry Smith - ovl - the amount of overlap between subdomains (ovl >= 0, default value = 1)
10224b9ad928SBarry Smith 
10234b9ad928SBarry Smith   Options Database Key:
10244b9ad928SBarry Smith . -pc_asm_overlap <ovl> - Sets overlap
10254b9ad928SBarry Smith 
102620f4b53cSBarry Smith   Level: intermediate
102720f4b53cSBarry Smith 
10284b9ad928SBarry Smith   Notes:
1029f1580f4eSBarry Smith   By default the `PCASM` preconditioner uses 1 block per processor.  To use
1030f1580f4eSBarry Smith   multiple blocks per perocessor, see `PCASMSetTotalSubdomains()` and
1031f1580f4eSBarry Smith   `PCASMSetLocalSubdomains()` (and the option -pc_asm_blocks <blks>).
10324b9ad928SBarry Smith 
10334b9ad928SBarry Smith   The overlap defaults to 1, so if one desires that no additional
10344b9ad928SBarry Smith   overlap be computed beyond what may have been set with a call to
1035f1580f4eSBarry Smith   `PCASMSetTotalSubdomains()` or `PCASMSetLocalSubdomains()`, then ovl
10364b9ad928SBarry Smith   must be set to be 0.  In particular, if one does not explicitly set
10374b9ad928SBarry Smith   the subdomains an application code, then all overlap would be computed
1038f1580f4eSBarry Smith   internally by PETSc, and using an overlap of 0 would result in an `PCASM`
10394b9ad928SBarry Smith   variant that is equivalent to the block Jacobi preconditioner.
10404b9ad928SBarry Smith 
1041f1ee410cSBarry Smith   The default algorithm used by PETSc to increase overlap is fast, but not scalable,
1042f1ee410cSBarry Smith   use the option -mat_increase_overlap_scalable when the problem and number of processes is large.
1043f1ee410cSBarry Smith 
10442fe279fdSBarry Smith   One can define initial index sets with any overlap via
1045f1580f4eSBarry Smith   `PCASMSetLocalSubdomains()`; the routine
1046f1580f4eSBarry Smith   `PCASMSetOverlap()` merely allows PETSc to extend that overlap further
10474b9ad928SBarry Smith   if desired.
10484b9ad928SBarry Smith 
1049562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMSetTotalSubdomains()`, `PCASMSetLocalSubdomains()`, `PCASMGetSubKSP()`,
1050f1580f4eSBarry Smith           `PCASMCreateSubdomains2D()`, `PCASMGetLocalSubdomains()`, `MatIncreaseOverlap()`, `PCGASM`
10514b9ad928SBarry Smith @*/
PCASMSetOverlap(PC pc,PetscInt ovl)1052d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMSetOverlap(PC pc, PetscInt ovl)
1053d71ae5a4SJacob Faibussowitsch {
10544b9ad928SBarry Smith   PetscFunctionBegin;
10550700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1056c5eb9154SBarry Smith   PetscValidLogicalCollectiveInt(pc, ovl, 2);
1057cac4c232SBarry Smith   PetscTryMethod(pc, "PCASMSetOverlap_C", (PC, PetscInt), (pc, ovl));
10583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10594b9ad928SBarry Smith }
10604b9ad928SBarry Smith 
10614b9ad928SBarry Smith /*@
10624b9ad928SBarry Smith   PCASMSetType - Sets the type of restriction and interpolation used
1063f1580f4eSBarry Smith   for local problems in the additive Schwarz method, `PCASM`.
10644b9ad928SBarry Smith 
1065c3339decSBarry Smith   Logically Collective
10664b9ad928SBarry Smith 
10674b9ad928SBarry Smith   Input Parameters:
10684b9ad928SBarry Smith + pc   - the preconditioner context
1069f1580f4eSBarry Smith - type - variant of `PCASM`, one of
10704b9ad928SBarry Smith .vb
10714b9ad928SBarry Smith       PC_ASM_BASIC       - full interpolation and restriction
107282b5ce2aSStefano Zampini       PC_ASM_RESTRICT    - full restriction, local processor interpolation (default)
10734b9ad928SBarry Smith       PC_ASM_INTERPOLATE - full interpolation, local processor restriction
10744b9ad928SBarry Smith       PC_ASM_NONE        - local processor restriction and interpolation
10754b9ad928SBarry Smith .ve
10764b9ad928SBarry Smith 
10774b9ad928SBarry Smith   Options Database Key:
1078f1580f4eSBarry Smith . -pc_asm_type [basic,restrict,interpolate,none] - Sets `PCASMType`
10794b9ad928SBarry Smith 
108020f4b53cSBarry Smith   Level: intermediate
108120f4b53cSBarry Smith 
1082f1580f4eSBarry Smith   Note:
1083f1580f4eSBarry Smith   if the is_local arguments are passed to `PCASMSetLocalSubdomains()` then they are used when `PC_ASM_RESTRICT` has been selected
1084f1ee410cSBarry Smith   to limit the local processor interpolation
1085f1ee410cSBarry Smith 
1086562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMSetTotalSubdomains()`, `PCASMGetSubKSP()`,
1087f1580f4eSBarry Smith           `PCASMCreateSubdomains2D()`, `PCASMType`, `PCASMSetLocalType()`, `PCASMGetLocalType()`, `PCGASM`
10884b9ad928SBarry Smith @*/
PCASMSetType(PC pc,PCASMType type)1089d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMSetType(PC pc, PCASMType type)
1090d71ae5a4SJacob Faibussowitsch {
10914b9ad928SBarry Smith   PetscFunctionBegin;
10920700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1093c5eb9154SBarry Smith   PetscValidLogicalCollectiveEnum(pc, type, 2);
1094cac4c232SBarry Smith   PetscTryMethod(pc, "PCASMSetType_C", (PC, PCASMType), (pc, type));
10953ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
10964b9ad928SBarry Smith }
10974b9ad928SBarry Smith 
1098c60c7ad4SBarry Smith /*@
1099c60c7ad4SBarry Smith   PCASMGetType - Gets the type of restriction and interpolation used
1100f1580f4eSBarry Smith   for local problems in the additive Schwarz method, `PCASM`.
1101c60c7ad4SBarry Smith 
1102c3339decSBarry Smith   Logically Collective
1103c60c7ad4SBarry Smith 
1104c60c7ad4SBarry Smith   Input Parameter:
1105c60c7ad4SBarry Smith . pc - the preconditioner context
1106c60c7ad4SBarry Smith 
1107c60c7ad4SBarry Smith   Output Parameter:
1108f1580f4eSBarry Smith . type - variant of `PCASM`, one of
1109c60c7ad4SBarry Smith .vb
1110c60c7ad4SBarry Smith       PC_ASM_BASIC       - full interpolation and restriction
1111c60c7ad4SBarry Smith       PC_ASM_RESTRICT    - full restriction, local processor interpolation
1112c60c7ad4SBarry Smith       PC_ASM_INTERPOLATE - full interpolation, local processor restriction
1113c60c7ad4SBarry Smith       PC_ASM_NONE        - local processor restriction and interpolation
1114c60c7ad4SBarry Smith .ve
1115c60c7ad4SBarry Smith 
1116c60c7ad4SBarry Smith   Options Database Key:
1117f1580f4eSBarry Smith . -pc_asm_type [basic,restrict,interpolate,none] - Sets `PCASM` type
1118c60c7ad4SBarry Smith 
1119c60c7ad4SBarry Smith   Level: intermediate
1120c60c7ad4SBarry Smith 
1121562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMSetTotalSubdomains()`, `PCASMGetSubKSP()`, `PCGASM`,
1122db781477SPatrick Sanan           `PCASMCreateSubdomains2D()`, `PCASMType`, `PCASMSetType()`, `PCASMSetLocalType()`, `PCASMGetLocalType()`
1123c60c7ad4SBarry Smith @*/
PCASMGetType(PC pc,PCASMType * type)1124d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMGetType(PC pc, PCASMType *type)
1125d71ae5a4SJacob Faibussowitsch {
1126c60c7ad4SBarry Smith   PetscFunctionBegin;
1127c60c7ad4SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1128cac4c232SBarry Smith   PetscUseMethod(pc, "PCASMGetType_C", (PC, PCASMType *), (pc, type));
11293ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1130c60c7ad4SBarry Smith }
1131c60c7ad4SBarry Smith 
113212cd4985SMatthew G. Knepley /*@
1133f1580f4eSBarry Smith   PCASMSetLocalType - Sets the type of composition used for local problems in the additive Schwarz method, `PCASM`.
113412cd4985SMatthew G. Knepley 
1135c3339decSBarry Smith   Logically Collective
113612cd4985SMatthew G. Knepley 
113712cd4985SMatthew G. Knepley   Input Parameters:
113812cd4985SMatthew G. Knepley + pc   - the preconditioner context
113912cd4985SMatthew G. Knepley - type - type of composition, one of
114012cd4985SMatthew G. Knepley .vb
114112cd4985SMatthew G. Knepley   PC_COMPOSITE_ADDITIVE       - local additive combination
114212cd4985SMatthew G. Knepley   PC_COMPOSITE_MULTIPLICATIVE - local multiplicative combination
114312cd4985SMatthew G. Knepley .ve
114412cd4985SMatthew G. Knepley 
114512cd4985SMatthew G. Knepley   Options Database Key:
114612cd4985SMatthew G. Knepley . -pc_asm_local_type [additive,multiplicative] - Sets local solver composition type
114712cd4985SMatthew G. Knepley 
114812cd4985SMatthew G. Knepley   Level: intermediate
114912cd4985SMatthew G. Knepley 
1150562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMSetType()`, `PCASMGetType()`, `PCASMGetLocalType()`, `PCASMType`, `PCCompositeType`
115112cd4985SMatthew G. Knepley @*/
PCASMSetLocalType(PC pc,PCCompositeType type)1152d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMSetLocalType(PC pc, PCCompositeType type)
1153d71ae5a4SJacob Faibussowitsch {
115412cd4985SMatthew G. Knepley   PetscFunctionBegin;
115512cd4985SMatthew G. Knepley   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
115612cd4985SMatthew G. Knepley   PetscValidLogicalCollectiveEnum(pc, type, 2);
1157cac4c232SBarry Smith   PetscTryMethod(pc, "PCASMSetLocalType_C", (PC, PCCompositeType), (pc, type));
11583ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
115912cd4985SMatthew G. Knepley }
116012cd4985SMatthew G. Knepley 
116112cd4985SMatthew G. Knepley /*@
1162f1580f4eSBarry Smith   PCASMGetLocalType - Gets the type of composition used for local problems in the additive Schwarz method, `PCASM`.
116312cd4985SMatthew G. Knepley 
1164c3339decSBarry Smith   Logically Collective
116512cd4985SMatthew G. Knepley 
116612cd4985SMatthew G. Knepley   Input Parameter:
116712cd4985SMatthew G. Knepley . pc - the preconditioner context
116812cd4985SMatthew G. Knepley 
116912cd4985SMatthew G. Knepley   Output Parameter:
117012cd4985SMatthew G. Knepley . type - type of composition, one of
117112cd4985SMatthew G. Knepley .vb
117212cd4985SMatthew G. Knepley   PC_COMPOSITE_ADDITIVE       - local additive combination
117312cd4985SMatthew G. Knepley   PC_COMPOSITE_MULTIPLICATIVE - local multiplicative combination
117412cd4985SMatthew G. Knepley .ve
117512cd4985SMatthew G. Knepley 
117612cd4985SMatthew G. Knepley   Options Database Key:
117712cd4985SMatthew G. Knepley . -pc_asm_local_type [additive,multiplicative] - Sets local solver composition type
117812cd4985SMatthew G. Knepley 
117912cd4985SMatthew G. Knepley   Level: intermediate
118012cd4985SMatthew G. Knepley 
11810241b274SPierre Jolivet .seealso: [](ch_ksp), `PCASM`, `PCASMSetType()`, `PCASMGetType()`, `PCASMSetLocalType()`, `PCASMType`, `PCCompositeType`
118212cd4985SMatthew G. Knepley @*/
PCASMGetLocalType(PC pc,PCCompositeType * type)1183d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMGetLocalType(PC pc, PCCompositeType *type)
1184d71ae5a4SJacob Faibussowitsch {
118512cd4985SMatthew G. Knepley   PetscFunctionBegin;
118612cd4985SMatthew G. Knepley   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
11874f572ea9SToby Isaac   PetscAssertPointer(type, 2);
1188cac4c232SBarry Smith   PetscUseMethod(pc, "PCASMGetLocalType_C", (PC, PCCompositeType *), (pc, type));
11893ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
119012cd4985SMatthew G. Knepley }
119112cd4985SMatthew G. Knepley 
11926ed231c7SMatthew Knepley /*@
11936ed231c7SMatthew Knepley   PCASMSetSortIndices - Determines whether subdomain indices are sorted.
11946ed231c7SMatthew Knepley 
1195c3339decSBarry Smith   Logically Collective
11966ed231c7SMatthew Knepley 
11976ed231c7SMatthew Knepley   Input Parameters:
11986ed231c7SMatthew Knepley + pc     - the preconditioner context
11996ed231c7SMatthew Knepley - doSort - sort the subdomain indices
12006ed231c7SMatthew Knepley 
12016ed231c7SMatthew Knepley   Level: intermediate
12026ed231c7SMatthew Knepley 
1203562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMSetLocalSubdomains()`, `PCASMSetTotalSubdomains()`, `PCASMGetSubKSP()`,
1204db781477SPatrick Sanan           `PCASMCreateSubdomains2D()`
12056ed231c7SMatthew Knepley @*/
PCASMSetSortIndices(PC pc,PetscBool doSort)1206d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMSetSortIndices(PC pc, PetscBool doSort)
1207d71ae5a4SJacob Faibussowitsch {
12086ed231c7SMatthew Knepley   PetscFunctionBegin;
12090700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1210acfcf0e5SJed Brown   PetscValidLogicalCollectiveBool(pc, doSort, 2);
1211cac4c232SBarry Smith   PetscTryMethod(pc, "PCASMSetSortIndices_C", (PC, PetscBool), (pc, doSort));
12123ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
12136ed231c7SMatthew Knepley }
12146ed231c7SMatthew Knepley 
12154b9ad928SBarry Smith /*@C
1216f1580f4eSBarry Smith   PCASMGetSubKSP - Gets the local `KSP` contexts for all blocks on
12174b9ad928SBarry Smith   this processor.
12184b9ad928SBarry Smith 
1219c3339decSBarry Smith   Collective iff first_local is requested
12204b9ad928SBarry Smith 
12214b9ad928SBarry Smith   Input Parameter:
12224b9ad928SBarry Smith . pc - the preconditioner context
12234b9ad928SBarry Smith 
12244b9ad928SBarry Smith   Output Parameters:
1225a3b724e8SBarry Smith + n_local     - the number of blocks on this processor or `NULL`
1226a3b724e8SBarry Smith . first_local - the global number of the first block on this processor or `NULL`, all processors must request or all must pass `NULL`
1227f1580f4eSBarry Smith - ksp         - the array of `KSP` contexts
12284b9ad928SBarry Smith 
122920f4b53cSBarry Smith   Level: advanced
123020f4b53cSBarry Smith 
1231f1580f4eSBarry Smith   Notes:
1232f1580f4eSBarry Smith   After `PCASMGetSubKSP()` the array of `KSP`s is not to be freed.
12334b9ad928SBarry Smith 
1234f1580f4eSBarry Smith   You must call `KSPSetUp()` before calling `PCASMGetSubKSP()`.
12354b9ad928SBarry Smith 
123636083efbSBarry Smith   Fortran Note:
123736083efbSBarry Smith   Call `PCASMRestoreSubKSP()` when access to the array of `KSP` is no longer needed
123836083efbSBarry Smith 
1239562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMSetTotalSubdomains()`, `PCASMSetOverlap()`,
1240db781477SPatrick Sanan           `PCASMCreateSubdomains2D()`,
12414b9ad928SBarry Smith @*/
PCASMGetSubKSP(PC pc,PetscInt * n_local,PetscInt * first_local,KSP * ksp[])1242d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMGetSubKSP(PC pc, PetscInt *n_local, PetscInt *first_local, KSP *ksp[])
1243d71ae5a4SJacob Faibussowitsch {
12444b9ad928SBarry Smith   PetscFunctionBegin;
12450700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1246cac4c232SBarry Smith   PetscUseMethod(pc, "PCASMGetSubKSP_C", (PC, PetscInt *, PetscInt *, KSP **), (pc, n_local, first_local, ksp));
12473ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
12484b9ad928SBarry Smith }
12494b9ad928SBarry Smith 
12504b9ad928SBarry Smith /*MC
12513b09bd56SBarry Smith    PCASM - Use the (restricted) additive Schwarz method, each block is (approximately) solved with
12521d27aa22SBarry Smith            its own `KSP` object, {cite}`dryja1987additive` and {cite}`1sbg`
12534b9ad928SBarry Smith 
12544b9ad928SBarry Smith    Options Database Keys:
12551d27aa22SBarry Smith +  -pc_asm_blocks <blks>                          - Sets total blocks. Defaults to one block per MPI process.
12564b9ad928SBarry Smith .  -pc_asm_overlap <ovl>                          - Sets overlap
1257f1580f4eSBarry Smith .  -pc_asm_type [basic,restrict,interpolate,none] - Sets `PCASMType`, default is restrict. See `PCASMSetType()`
125873ff1848SBarry Smith .  -pc_asm_dm_subdomains <bool>                   - use subdomains defined by the `DM` with `DMCreateDomainDecomposition()`
1259f1580f4eSBarry Smith -  -pc_asm_local_type [additive, multiplicative]  - Sets `PCCompositeType`, default is additive. See `PCASMSetLocalType()`
12604b9ad928SBarry Smith 
12614b9ad928SBarry Smith    Level: beginner
12624b9ad928SBarry Smith 
1263f1580f4eSBarry Smith    Notes:
1264f1580f4eSBarry Smith    If you run with, for example, 3 blocks on 1 processor or 3 blocks on 3 processors you
126573ff1848SBarry Smith    will get a different convergence rate due to the default option of `-pc_asm_type restrict`. Use
126673ff1848SBarry Smith    `-pc_asm_type basic` to get the same convergence behavior
1267f1580f4eSBarry Smith 
1268f1580f4eSBarry Smith    Each processor can have one or more blocks, but a block cannot be shared by more
1269f1580f4eSBarry Smith    than one processor. Use `PCGASM` for subdomains shared by multiple processes.
1270f1580f4eSBarry Smith 
127173ff1848SBarry Smith    To set options on the solvers for each block append `-sub_` to all the `KSP`, and `PC`
127273ff1848SBarry Smith    options database keys. For example, `-sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly`
1273f1580f4eSBarry Smith 
1274f1580f4eSBarry Smith    To set the options on the solvers separate for each block call `PCASMGetSubKSP()`
1275f1580f4eSBarry Smith    and set the options directly on the resulting `KSP` object (you can access its `PC` with `KSPGetPC()`)
1276f1580f4eSBarry Smith 
127773ff1848SBarry Smith    If the `PC` has an associated `DM`, then, by default, `DMCreateDomainDecomposition()` is used to create the subdomains
127873ff1848SBarry Smith 
1279562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCASMType`, `PCCompositeType`,
1280db781477SPatrick Sanan           `PCBJACOBI`, `PCASMGetSubKSP()`, `PCASMSetLocalSubdomains()`, `PCASMType`, `PCASMGetType()`, `PCASMSetLocalType()`, `PCASMGetLocalType()`
1281db781477SPatrick Sanan           `PCASMSetTotalSubdomains()`, `PCSetModifySubMatrices()`, `PCASMSetOverlap()`, `PCASMSetType()`, `PCCompositeType`
12824b9ad928SBarry Smith M*/
12834b9ad928SBarry Smith 
PCCreate_ASM(PC pc)1284d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_ASM(PC pc)
1285d71ae5a4SJacob Faibussowitsch {
12864b9ad928SBarry Smith   PC_ASM *osm;
12874b9ad928SBarry Smith 
12884b9ad928SBarry Smith   PetscFunctionBegin;
12894dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&osm));
12902fa5cd67SKarl Rupp 
12914b9ad928SBarry Smith   osm->n             = PETSC_DECIDE;
12924b9ad928SBarry Smith   osm->n_local       = 0;
12932b837212SDmitry Karpeev   osm->n_local_true  = PETSC_DECIDE;
12944b9ad928SBarry Smith   osm->overlap       = 1;
12950a545947SLisandro Dalcin   osm->ksp           = NULL;
12960a545947SLisandro Dalcin   osm->restriction   = NULL;
12970a545947SLisandro Dalcin   osm->lprolongation = NULL;
12980a545947SLisandro Dalcin   osm->lrestriction  = NULL;
12990a545947SLisandro Dalcin   osm->x             = NULL;
13000a545947SLisandro Dalcin   osm->y             = NULL;
13010a545947SLisandro Dalcin   osm->is            = NULL;
13020a545947SLisandro Dalcin   osm->is_local      = NULL;
13030a545947SLisandro Dalcin   osm->mat           = NULL;
13040a545947SLisandro Dalcin   osm->pmat          = NULL;
13054b9ad928SBarry Smith   osm->type          = PC_ASM_RESTRICT;
130612cd4985SMatthew G. Knepley   osm->loctype       = PC_COMPOSITE_ADDITIVE;
13076ed231c7SMatthew Knepley   osm->sort_indices  = PETSC_TRUE;
1308d709ea83SDmitry Karpeev   osm->dm_subdomains = PETSC_FALSE;
130980ec0b7dSPatrick Sanan   osm->sub_mat_type  = NULL;
13104b9ad928SBarry Smith 
131192bb6962SLisandro Dalcin   pc->data                   = (void *)osm;
13124b9ad928SBarry Smith   pc->ops->apply             = PCApply_ASM;
131348e38a8aSPierre Jolivet   pc->ops->matapply          = PCMatApply_ASM;
13144b9ad928SBarry Smith   pc->ops->applytranspose    = PCApplyTranspose_ASM;
13154dbf25a8SPierre Jolivet   pc->ops->matapplytranspose = PCMatApplyTranspose_ASM;
13164b9ad928SBarry Smith   pc->ops->setup             = PCSetUp_ASM;
1317e91c6855SBarry Smith   pc->ops->reset             = PCReset_ASM;
13184b9ad928SBarry Smith   pc->ops->destroy           = PCDestroy_ASM;
13194b9ad928SBarry Smith   pc->ops->setfromoptions    = PCSetFromOptions_ASM;
13204b9ad928SBarry Smith   pc->ops->setuponblocks     = PCSetUpOnBlocks_ASM;
13214b9ad928SBarry Smith   pc->ops->view              = PCView_ASM;
13220a545947SLisandro Dalcin   pc->ops->applyrichardson   = NULL;
13234b9ad928SBarry Smith 
13249566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetLocalSubdomains_C", PCASMSetLocalSubdomains_ASM));
13259566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetTotalSubdomains_C", PCASMSetTotalSubdomains_ASM));
13269566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetOverlap_C", PCASMSetOverlap_ASM));
13279566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetType_C", PCASMSetType_ASM));
13289566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMGetType_C", PCASMGetType_ASM));
13299566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetLocalType_C", PCASMSetLocalType_ASM));
13309566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMGetLocalType_C", PCASMGetLocalType_ASM));
13319566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetSortIndices_C", PCASMSetSortIndices_ASM));
13329566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMGetSubKSP_C", PCASMGetSubKSP_ASM));
13339566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMGetSubMatType_C", PCASMGetSubMatType_ASM));
13349566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetSubMatType_C", PCASMSetSubMatType_ASM));
13353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13364b9ad928SBarry Smith }
13374b9ad928SBarry Smith 
133892bb6962SLisandro Dalcin /*@C
133992bb6962SLisandro Dalcin   PCASMCreateSubdomains - Creates the index sets for the overlapping Schwarz
1340f1580f4eSBarry Smith   preconditioner, `PCASM`,  for any problem on a general grid.
134192bb6962SLisandro Dalcin 
134292bb6962SLisandro Dalcin   Collective
134392bb6962SLisandro Dalcin 
134492bb6962SLisandro Dalcin   Input Parameters:
134592bb6962SLisandro Dalcin + A - The global matrix operator
134692bb6962SLisandro Dalcin - n - the number of local blocks
134792bb6962SLisandro Dalcin 
13482fe279fdSBarry Smith   Output Parameter:
134992bb6962SLisandro Dalcin . outis - the array of index sets defining the subdomains
135092bb6962SLisandro Dalcin 
135192bb6962SLisandro Dalcin   Level: advanced
135292bb6962SLisandro Dalcin 
1353f1580f4eSBarry Smith   Note:
1354f1580f4eSBarry Smith   This generates nonoverlapping subdomains; the `PCASM` will generate the overlap
1355f1580f4eSBarry Smith   from these if you use `PCASMSetLocalSubdomains()`
13567d6bfa3bSBarry Smith 
1357562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMSetLocalSubdomains()`, `PCASMDestroySubdomains()`
135892bb6962SLisandro Dalcin @*/
PCASMCreateSubdomains(Mat A,PetscInt n,IS * outis[])1359d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMCreateSubdomains(Mat A, PetscInt n, IS *outis[])
1360d71ae5a4SJacob Faibussowitsch {
136192bb6962SLisandro Dalcin   MatPartitioning mpart;
136292bb6962SLisandro Dalcin   const char     *prefix;
136392bb6962SLisandro Dalcin   PetscInt        i, j, rstart, rend, bs;
1364976e8c5aSLisandro Dalcin   PetscBool       hasop, isbaij = PETSC_FALSE, foundpart = PETSC_FALSE;
13650298fd71SBarry Smith   Mat             Ad = NULL, adj;
136692bb6962SLisandro Dalcin   IS              ispart, isnumb, *is;
136792bb6962SLisandro Dalcin 
136892bb6962SLisandro Dalcin   PetscFunctionBegin;
13690700a824SBarry Smith   PetscValidHeaderSpecific(A, MAT_CLASSID, 1);
13704f572ea9SToby Isaac   PetscAssertPointer(outis, 3);
137163a3b9bcSJacob Faibussowitsch   PetscCheck(n >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "number of local blocks must be > 0, n = %" PetscInt_FMT, n);
137292bb6962SLisandro Dalcin 
137392bb6962SLisandro Dalcin   /* Get prefix, row distribution, and block size */
13749566063dSJacob Faibussowitsch   PetscCall(MatGetOptionsPrefix(A, &prefix));
13759566063dSJacob Faibussowitsch   PetscCall(MatGetOwnershipRange(A, &rstart, &rend));
13769566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(A, &bs));
137763a3b9bcSJacob 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);
137865e19b50SBarry Smith 
137992bb6962SLisandro Dalcin   /* Get diagonal block from matrix if possible */
13809566063dSJacob Faibussowitsch   PetscCall(MatHasOperation(A, MATOP_GET_DIAGONAL_BLOCK, &hasop));
138148a46eb9SPierre Jolivet   if (hasop) PetscCall(MatGetDiagonalBlock(A, &Ad));
138292bb6962SLisandro Dalcin   if (Ad) {
13839566063dSJacob Faibussowitsch     PetscCall(PetscObjectBaseTypeCompare((PetscObject)Ad, MATSEQBAIJ, &isbaij));
13849566063dSJacob Faibussowitsch     if (!isbaij) PetscCall(PetscObjectBaseTypeCompare((PetscObject)Ad, MATSEQSBAIJ, &isbaij));
138592bb6962SLisandro Dalcin   }
138692bb6962SLisandro Dalcin   if (Ad && n > 1) {
1387ace3abfcSBarry Smith     PetscBool match, done;
138892bb6962SLisandro Dalcin     /* Try to setup a good matrix partitioning if available */
13899566063dSJacob Faibussowitsch     PetscCall(MatPartitioningCreate(PETSC_COMM_SELF, &mpart));
13909566063dSJacob Faibussowitsch     PetscCall(PetscObjectSetOptionsPrefix((PetscObject)mpart, prefix));
13919566063dSJacob Faibussowitsch     PetscCall(MatPartitioningSetFromOptions(mpart));
13929566063dSJacob Faibussowitsch     PetscCall(PetscObjectTypeCompare((PetscObject)mpart, MATPARTITIONINGCURRENT, &match));
139348a46eb9SPierre Jolivet     if (!match) PetscCall(PetscObjectTypeCompare((PetscObject)mpart, MATPARTITIONINGSQUARE, &match));
139492bb6962SLisandro Dalcin     if (!match) { /* assume a "good" partitioner is available */
13951a83f524SJed Brown       PetscInt        na;
13961a83f524SJed Brown       const PetscInt *ia, *ja;
13979566063dSJacob Faibussowitsch       PetscCall(MatGetRowIJ(Ad, 0, PETSC_TRUE, isbaij, &na, &ia, &ja, &done));
139892bb6962SLisandro Dalcin       if (done) {
139992bb6962SLisandro Dalcin         /* Build adjacency matrix by hand. Unfortunately a call to
140092bb6962SLisandro Dalcin            MatConvert(Ad,MATMPIADJ,MAT_INITIAL_MATRIX,&adj) will
140192bb6962SLisandro Dalcin            remove the block-aij structure and we cannot expect
140292bb6962SLisandro Dalcin            MatPartitioning to split vertices as we need */
14030a545947SLisandro Dalcin         PetscInt        i, j, len, nnz, cnt, *iia = NULL, *jja = NULL;
14041a83f524SJed Brown         const PetscInt *row;
140592bb6962SLisandro Dalcin         nnz = 0;
140692bb6962SLisandro Dalcin         for (i = 0; i < na; i++) { /* count number of nonzeros */
140792bb6962SLisandro Dalcin           len = ia[i + 1] - ia[i];
140892bb6962SLisandro Dalcin           row = ja + ia[i];
140992bb6962SLisandro Dalcin           for (j = 0; j < len; j++) {
141092bb6962SLisandro Dalcin             if (row[j] == i) { /* don't count diagonal */
14119371c9d4SSatish Balay               len--;
14129371c9d4SSatish Balay               break;
141392bb6962SLisandro Dalcin             }
141492bb6962SLisandro Dalcin           }
141592bb6962SLisandro Dalcin           nnz += len;
141692bb6962SLisandro Dalcin         }
14179566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(na + 1, &iia));
14189566063dSJacob Faibussowitsch         PetscCall(PetscMalloc1(nnz, &jja));
141992bb6962SLisandro Dalcin         nnz    = 0;
142092bb6962SLisandro Dalcin         iia[0] = 0;
142192bb6962SLisandro Dalcin         for (i = 0; i < na; i++) { /* fill adjacency */
142292bb6962SLisandro Dalcin           cnt = 0;
142392bb6962SLisandro Dalcin           len = ia[i + 1] - ia[i];
142492bb6962SLisandro Dalcin           row = ja + ia[i];
142592bb6962SLisandro Dalcin           for (j = 0; j < len; j++) {
142692bb6962SLisandro Dalcin             if (row[j] != i) { /* if not diagonal */
142792bb6962SLisandro Dalcin               jja[nnz + cnt++] = row[j];
142892bb6962SLisandro Dalcin             }
142992bb6962SLisandro Dalcin           }
143092bb6962SLisandro Dalcin           nnz += cnt;
143192bb6962SLisandro Dalcin           iia[i + 1] = nnz;
143292bb6962SLisandro Dalcin         }
143392bb6962SLisandro Dalcin         /* Partitioning of the adjacency matrix */
14349566063dSJacob Faibussowitsch         PetscCall(MatCreateMPIAdj(PETSC_COMM_SELF, na, na, iia, jja, NULL, &adj));
14359566063dSJacob Faibussowitsch         PetscCall(MatPartitioningSetAdjacency(mpart, adj));
14369566063dSJacob Faibussowitsch         PetscCall(MatPartitioningSetNParts(mpart, n));
14379566063dSJacob Faibussowitsch         PetscCall(MatPartitioningApply(mpart, &ispart));
14389566063dSJacob Faibussowitsch         PetscCall(ISPartitioningToNumbering(ispart, &isnumb));
14399566063dSJacob Faibussowitsch         PetscCall(MatDestroy(&adj));
144092bb6962SLisandro Dalcin         foundpart = PETSC_TRUE;
144192bb6962SLisandro Dalcin       }
14429566063dSJacob Faibussowitsch       PetscCall(MatRestoreRowIJ(Ad, 0, PETSC_TRUE, isbaij, &na, &ia, &ja, &done));
144392bb6962SLisandro Dalcin     }
14449566063dSJacob Faibussowitsch     PetscCall(MatPartitioningDestroy(&mpart));
144592bb6962SLisandro Dalcin   }
144692bb6962SLisandro Dalcin 
14479566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n, &is));
144892bb6962SLisandro Dalcin   *outis = is;
144992bb6962SLisandro Dalcin 
145092bb6962SLisandro Dalcin   if (!foundpart) {
145192bb6962SLisandro Dalcin     /* Partitioning by contiguous chunks of rows */
145292bb6962SLisandro Dalcin 
145392bb6962SLisandro Dalcin     PetscInt mbs   = (rend - rstart) / bs;
145492bb6962SLisandro Dalcin     PetscInt start = rstart;
145592bb6962SLisandro Dalcin     for (i = 0; i < n; i++) {
145692bb6962SLisandro Dalcin       PetscInt count = (mbs / n + ((mbs % n) > i)) * bs;
14579566063dSJacob Faibussowitsch       PetscCall(ISCreateStride(PETSC_COMM_SELF, count, start, 1, &is[i]));
145892bb6962SLisandro Dalcin       start += count;
145992bb6962SLisandro Dalcin     }
146092bb6962SLisandro Dalcin 
146192bb6962SLisandro Dalcin   } else {
146292bb6962SLisandro Dalcin     /* Partitioning by adjacency of diagonal block  */
146392bb6962SLisandro Dalcin 
146492bb6962SLisandro Dalcin     const PetscInt *numbering;
146592bb6962SLisandro Dalcin     PetscInt       *count, nidx, *indices, *newidx, start = 0;
146692bb6962SLisandro Dalcin     /* Get node count in each partition */
14679566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(n, &count));
14689566063dSJacob Faibussowitsch     PetscCall(ISPartitioningCount(ispart, n, count));
146992bb6962SLisandro Dalcin     if (isbaij && bs > 1) { /* adjust for the block-aij case */
147092bb6962SLisandro Dalcin       for (i = 0; i < n; i++) count[i] *= bs;
147192bb6962SLisandro Dalcin     }
147292bb6962SLisandro Dalcin     /* Build indices from node numbering */
14739566063dSJacob Faibussowitsch     PetscCall(ISGetLocalSize(isnumb, &nidx));
14749566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(nidx, &indices));
147592bb6962SLisandro Dalcin     for (i = 0; i < nidx; i++) indices[i] = i; /* needs to be initialized */
14769566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(isnumb, &numbering));
14779566063dSJacob Faibussowitsch     PetscCall(PetscSortIntWithPermutation(nidx, numbering, indices));
14789566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(isnumb, &numbering));
147992bb6962SLisandro Dalcin     if (isbaij && bs > 1) { /* adjust for the block-aij case */
14809566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(nidx * bs, &newidx));
14812fa5cd67SKarl Rupp       for (i = 0; i < nidx; i++) {
14822fa5cd67SKarl Rupp         for (j = 0; j < bs; j++) newidx[i * bs + j] = indices[i] * bs + j;
14832fa5cd67SKarl Rupp       }
14849566063dSJacob Faibussowitsch       PetscCall(PetscFree(indices));
148592bb6962SLisandro Dalcin       nidx *= bs;
148692bb6962SLisandro Dalcin       indices = newidx;
148792bb6962SLisandro Dalcin     }
148892bb6962SLisandro Dalcin     /* Shift to get global indices */
148992bb6962SLisandro Dalcin     for (i = 0; i < nidx; i++) indices[i] += rstart;
149092bb6962SLisandro Dalcin 
149192bb6962SLisandro Dalcin     /* Build the index sets for each block */
149292bb6962SLisandro Dalcin     for (i = 0; i < n; i++) {
14939566063dSJacob Faibussowitsch       PetscCall(ISCreateGeneral(PETSC_COMM_SELF, count[i], &indices[start], PETSC_COPY_VALUES, &is[i]));
14949566063dSJacob Faibussowitsch       PetscCall(ISSort(is[i]));
149592bb6962SLisandro Dalcin       start += count[i];
149692bb6962SLisandro Dalcin     }
149792bb6962SLisandro Dalcin 
14989566063dSJacob Faibussowitsch     PetscCall(PetscFree(count));
14999566063dSJacob Faibussowitsch     PetscCall(PetscFree(indices));
15009566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&isnumb));
15019566063dSJacob Faibussowitsch     PetscCall(ISDestroy(&ispart));
150292bb6962SLisandro Dalcin   }
15033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
150492bb6962SLisandro Dalcin }
150592bb6962SLisandro Dalcin 
150692bb6962SLisandro Dalcin /*@C
150792bb6962SLisandro Dalcin   PCASMDestroySubdomains - Destroys the index sets created with
1508f1580f4eSBarry Smith   `PCASMCreateSubdomains()`. Should be called after setting subdomains with `PCASMSetLocalSubdomains()`.
150992bb6962SLisandro Dalcin 
151092bb6962SLisandro Dalcin   Collective
151192bb6962SLisandro Dalcin 
151292bb6962SLisandro Dalcin   Input Parameters:
151392bb6962SLisandro Dalcin + n        - the number of index sets
15142b691e39SMatthew Knepley . is       - the array of index sets
15152fe279fdSBarry Smith - is_local - the array of local index sets, can be `NULL`
151692bb6962SLisandro Dalcin 
151792bb6962SLisandro Dalcin   Level: advanced
151892bb6962SLisandro Dalcin 
1519ce78bad3SBarry Smith   Developer Note:
1520ce78bad3SBarry Smith   The `IS` arguments should be a *[]
1521ce78bad3SBarry Smith 
1522562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMCreateSubdomains()`, `PCASMSetLocalSubdomains()`
152392bb6962SLisandro Dalcin @*/
PCASMDestroySubdomains(PetscInt n,IS * is[],IS * is_local[])1524ce78bad3SBarry Smith PetscErrorCode PCASMDestroySubdomains(PetscInt n, IS *is[], IS *is_local[])
1525d71ae5a4SJacob Faibussowitsch {
152692bb6962SLisandro Dalcin   PetscInt i;
15275fd66863SKarl Rupp 
152892bb6962SLisandro Dalcin   PetscFunctionBegin;
15293ba16761SJacob Faibussowitsch   if (n <= 0) PetscFunctionReturn(PETSC_SUCCESS);
1530ce78bad3SBarry Smith   if (*is) {
1531ce78bad3SBarry Smith     PetscAssertPointer(*is, 2);
1532ce78bad3SBarry Smith     for (i = 0; i < n; i++) PetscCall(ISDestroy(&(*is)[i]));
1533ce78bad3SBarry Smith     PetscCall(PetscFree(*is));
1534a35b7b57SDmitry Karpeev   }
1535ce78bad3SBarry Smith   if (is_local && *is_local) {
1536ce78bad3SBarry Smith     PetscAssertPointer(*is_local, 3);
1537ce78bad3SBarry Smith     for (i = 0; i < n; i++) PetscCall(ISDestroy(&(*is_local)[i]));
1538ce78bad3SBarry Smith     PetscCall(PetscFree(*is_local));
15392b691e39SMatthew Knepley   }
15403ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
154192bb6962SLisandro Dalcin }
154292bb6962SLisandro Dalcin 
15436141accfSBarry Smith /*@C
15444b9ad928SBarry Smith   PCASMCreateSubdomains2D - Creates the index sets for the overlapping Schwarz
1545f1580f4eSBarry Smith   preconditioner, `PCASM`, for a two-dimensional problem on a regular grid.
15464b9ad928SBarry Smith 
15474b9ad928SBarry Smith   Not Collective
15484b9ad928SBarry Smith 
15494b9ad928SBarry Smith   Input Parameters:
15506b867d5aSJose E. Roman + m       - the number of mesh points in the x direction
15516b867d5aSJose E. Roman . n       - the number of mesh points in the y direction
15526b867d5aSJose E. Roman . M       - the number of subdomains in the x direction
15536b867d5aSJose E. Roman . N       - the number of subdomains in the y direction
15544b9ad928SBarry Smith . dof     - degrees of freedom per node
15554b9ad928SBarry Smith - overlap - overlap in mesh lines
15564b9ad928SBarry Smith 
15574b9ad928SBarry Smith   Output Parameters:
15584b9ad928SBarry Smith + Nsub     - the number of subdomains created
15593d03bb48SJed Brown . is       - array of index sets defining overlapping (if overlap > 0) subdomains
15603d03bb48SJed Brown - is_local - array of index sets defining non-overlapping subdomains
15614b9ad928SBarry Smith 
156220f4b53cSBarry Smith   Level: advanced
156320f4b53cSBarry Smith 
15644b9ad928SBarry Smith   Note:
1565f1580f4eSBarry Smith   Presently `PCAMSCreateSubdomains2d()` is valid only for sequential
15664b9ad928SBarry Smith   preconditioners.  More general related routines are
1567f1580f4eSBarry Smith   `PCASMSetTotalSubdomains()` and `PCASMSetLocalSubdomains()`.
15684b9ad928SBarry Smith 
1569562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMSetTotalSubdomains()`, `PCASMSetLocalSubdomains()`, `PCASMGetSubKSP()`,
1570db781477SPatrick Sanan           `PCASMSetOverlap()`
15714b9ad928SBarry Smith @*/
PCASMCreateSubdomains2D(PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt dof,PetscInt overlap,PetscInt * Nsub,IS * is[],IS * is_local[])1572ce78bad3SBarry Smith PetscErrorCode PCASMCreateSubdomains2D(PetscInt m, PetscInt n, PetscInt M, PetscInt N, PetscInt dof, PetscInt overlap, PetscInt *Nsub, IS *is[], IS *is_local[])
1573d71ae5a4SJacob Faibussowitsch {
15743d03bb48SJed Brown   PetscInt i, j, height, width, ystart, xstart, yleft, yright, xleft, xright, loc_outer;
157513f74950SBarry Smith   PetscInt nidx, *idx, loc, ii, jj, count;
15764b9ad928SBarry Smith 
15774b9ad928SBarry Smith   PetscFunctionBegin;
15787827d75bSBarry Smith   PetscCheck(dof == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "dof must be 1");
15794b9ad928SBarry Smith 
15804b9ad928SBarry Smith   *Nsub = N * M;
15819566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(*Nsub, is));
15829566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(*Nsub, is_local));
15834b9ad928SBarry Smith   ystart    = 0;
15843d03bb48SJed Brown   loc_outer = 0;
15854b9ad928SBarry Smith   for (i = 0; i < N; i++) {
15864b9ad928SBarry Smith     height = n / N + ((n % N) > i); /* height of subdomain */
15877827d75bSBarry Smith     PetscCheck(height >= 2, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many N subdomains for mesh dimension n");
15889371c9d4SSatish Balay     yleft = ystart - overlap;
15899371c9d4SSatish Balay     if (yleft < 0) yleft = 0;
15909371c9d4SSatish Balay     yright = ystart + height + overlap;
15919371c9d4SSatish Balay     if (yright > n) yright = n;
15924b9ad928SBarry Smith     xstart = 0;
15934b9ad928SBarry Smith     for (j = 0; j < M; j++) {
15944b9ad928SBarry Smith       width = m / M + ((m % M) > j); /* width of subdomain */
15957827d75bSBarry Smith       PetscCheck(width >= 2, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many M subdomains for mesh dimension m");
15969371c9d4SSatish Balay       xleft = xstart - overlap;
15979371c9d4SSatish Balay       if (xleft < 0) xleft = 0;
15989371c9d4SSatish Balay       xright = xstart + width + overlap;
15999371c9d4SSatish Balay       if (xright > m) xright = m;
16004b9ad928SBarry Smith       nidx = (xright - xleft) * (yright - yleft);
16019566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(nidx, &idx));
16024b9ad928SBarry Smith       loc = 0;
16034b9ad928SBarry Smith       for (ii = yleft; ii < yright; ii++) {
16044b9ad928SBarry Smith         count = m * ii + xleft;
16052fa5cd67SKarl Rupp         for (jj = xleft; jj < xright; jj++) idx[loc++] = count++;
16064b9ad928SBarry Smith       }
16079566063dSJacob Faibussowitsch       PetscCall(ISCreateGeneral(PETSC_COMM_SELF, nidx, idx, PETSC_COPY_VALUES, (*is) + loc_outer));
16083d03bb48SJed Brown       if (overlap == 0) {
16099566063dSJacob Faibussowitsch         PetscCall(PetscObjectReference((PetscObject)(*is)[loc_outer]));
16102fa5cd67SKarl Rupp 
16113d03bb48SJed Brown         (*is_local)[loc_outer] = (*is)[loc_outer];
16123d03bb48SJed Brown       } else {
16133d03bb48SJed Brown         for (loc = 0, ii = ystart; ii < ystart + height; ii++) {
1614ad540459SPierre Jolivet           for (jj = xstart; jj < xstart + width; jj++) idx[loc++] = m * ii + jj;
16153d03bb48SJed Brown         }
16169566063dSJacob Faibussowitsch         PetscCall(ISCreateGeneral(PETSC_COMM_SELF, loc, idx, PETSC_COPY_VALUES, *is_local + loc_outer));
16173d03bb48SJed Brown       }
16189566063dSJacob Faibussowitsch       PetscCall(PetscFree(idx));
16194b9ad928SBarry Smith       xstart += width;
16203d03bb48SJed Brown       loc_outer++;
16214b9ad928SBarry Smith     }
16224b9ad928SBarry Smith     ystart += height;
16234b9ad928SBarry Smith   }
16249566063dSJacob Faibussowitsch   for (i = 0; i < *Nsub; i++) PetscCall(ISSort((*is)[i]));
16253ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
16264b9ad928SBarry Smith }
16274b9ad928SBarry Smith 
16284b9ad928SBarry Smith /*@C
16294b9ad928SBarry Smith   PCASMGetLocalSubdomains - Gets the local subdomains (for this processor
1630f1580f4eSBarry Smith   only) for the additive Schwarz preconditioner, `PCASM`.
16314b9ad928SBarry Smith 
1632ad4df100SBarry Smith   Not Collective
16334b9ad928SBarry Smith 
16344b9ad928SBarry Smith   Input Parameter:
16354b9ad928SBarry Smith . pc - the preconditioner context
16364b9ad928SBarry Smith 
16374b9ad928SBarry Smith   Output Parameters:
1638f17a6784SPierre Jolivet + n        - if requested, the number of subdomains for this processor (default value = 1)
1639f17a6784SPierre Jolivet . is       - if requested, the index sets that define the subdomains for this processor
164020f4b53cSBarry Smith - is_local - if requested, the index sets that define the local part of the subdomains for this processor (can be `NULL`)
164120f4b53cSBarry Smith 
164220f4b53cSBarry Smith   Level: advanced
16434b9ad928SBarry Smith 
1644f1580f4eSBarry Smith   Note:
1645f1580f4eSBarry Smith   The `IS` numbering is in the parallel, global numbering of the vector.
16464b9ad928SBarry Smith 
1647562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMSetTotalSubdomains()`, `PCASMSetOverlap()`, `PCASMGetSubKSP()`,
1648db781477SPatrick Sanan           `PCASMCreateSubdomains2D()`, `PCASMSetLocalSubdomains()`, `PCASMGetLocalSubmatrices()`
16494b9ad928SBarry Smith @*/
PCASMGetLocalSubdomains(PC pc,PetscInt * n,IS * is[],IS * is_local[])1650d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMGetLocalSubdomains(PC pc, PetscInt *n, IS *is[], IS *is_local[])
1651d71ae5a4SJacob Faibussowitsch {
16522a808120SBarry Smith   PC_ASM   *osm = (PC_ASM *)pc->data;
1653ace3abfcSBarry Smith   PetscBool match;
16544b9ad928SBarry Smith 
16554b9ad928SBarry Smith   PetscFunctionBegin;
16560700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
16574f572ea9SToby Isaac   if (n) PetscAssertPointer(n, 2);
16584f572ea9SToby Isaac   if (is) PetscAssertPointer(is, 3);
16594f572ea9SToby Isaac   if (is_local) PetscAssertPointer(is_local, 4);
16609566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCASM, &match));
166128b400f6SJacob Faibussowitsch   PetscCheck(match, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "PC is not a PCASM");
16624b9ad928SBarry Smith   if (n) *n = osm->n_local_true;
16634b9ad928SBarry Smith   if (is) *is = osm->is;
16642b691e39SMatthew Knepley   if (is_local) *is_local = osm->is_local;
16653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
16664b9ad928SBarry Smith }
16674b9ad928SBarry Smith 
16684b9ad928SBarry Smith /*@C
16694b9ad928SBarry Smith   PCASMGetLocalSubmatrices - Gets the local submatrices (for this processor
1670f1580f4eSBarry Smith   only) for the additive Schwarz preconditioner, `PCASM`.
16714b9ad928SBarry Smith 
1672ad4df100SBarry Smith   Not Collective
16734b9ad928SBarry Smith 
16744b9ad928SBarry Smith   Input Parameter:
16754b9ad928SBarry Smith . pc - the preconditioner context
16764b9ad928SBarry Smith 
16774b9ad928SBarry Smith   Output Parameters:
1678f17a6784SPierre Jolivet + n   - if requested, the number of matrices for this processor (default value = 1)
1679f17a6784SPierre Jolivet - mat - if requested, the matrices
16804b9ad928SBarry Smith 
16814b9ad928SBarry Smith   Level: advanced
16824b9ad928SBarry Smith 
168395452b02SPatrick Sanan   Notes:
1684f1580f4eSBarry Smith   Call after `PCSetUp()` (or `KSPSetUp()`) but before `PCApply()` and before `PCSetUpOnBlocks()`)
1685cf739d55SBarry Smith 
1686f1580f4eSBarry Smith   Usually one would use `PCSetModifySubMatrices()` to change the submatrices in building the preconditioner.
1687cf739d55SBarry Smith 
1688562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMSetTotalSubdomains()`, `PCASMSetOverlap()`, `PCASMGetSubKSP()`,
1689db781477SPatrick Sanan           `PCASMCreateSubdomains2D()`, `PCASMSetLocalSubdomains()`, `PCASMGetLocalSubdomains()`, `PCSetModifySubMatrices()`
16904b9ad928SBarry Smith @*/
PCASMGetLocalSubmatrices(PC pc,PetscInt * n,Mat * mat[])1691d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMGetLocalSubmatrices(PC pc, PetscInt *n, Mat *mat[])
1692d71ae5a4SJacob Faibussowitsch {
16934b9ad928SBarry Smith   PC_ASM   *osm;
1694ace3abfcSBarry Smith   PetscBool match;
16954b9ad928SBarry Smith 
16964b9ad928SBarry Smith   PetscFunctionBegin;
16970700a824SBarry Smith   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
16984f572ea9SToby Isaac   if (n) PetscAssertPointer(n, 2);
16994f572ea9SToby Isaac   if (mat) PetscAssertPointer(mat, 3);
170028b400f6SJacob Faibussowitsch   PetscCheck(pc->setupcalled, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must call after KSPSetUp() or PCSetUp().");
17019566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCASM, &match));
170292bb6962SLisandro Dalcin   if (!match) {
170392bb6962SLisandro Dalcin     if (n) *n = 0;
17040298fd71SBarry Smith     if (mat) *mat = NULL;
170592bb6962SLisandro Dalcin   } else {
17064b9ad928SBarry Smith     osm = (PC_ASM *)pc->data;
17074b9ad928SBarry Smith     if (n) *n = osm->n_local_true;
17084b9ad928SBarry Smith     if (mat) *mat = osm->pmat;
170992bb6962SLisandro Dalcin   }
17103ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
17114b9ad928SBarry Smith }
1712d709ea83SDmitry Karpeev 
1713d709ea83SDmitry Karpeev /*@
1714f1580f4eSBarry Smith   PCASMSetDMSubdomains - Indicates whether to use `DMCreateDomainDecomposition()` to define the subdomains, whenever possible.
1715f1ee410cSBarry Smith 
1716d709ea83SDmitry Karpeev   Logically Collective
1717d709ea83SDmitry Karpeev 
1718d8d19677SJose E. Roman   Input Parameters:
1719d709ea83SDmitry Karpeev + pc  - the preconditioner
1720f1580f4eSBarry Smith - flg - boolean indicating whether to use subdomains defined by the `DM`
1721d709ea83SDmitry Karpeev 
1722d709ea83SDmitry Karpeev   Options Database Key:
172373ff1848SBarry Smith . -pc_asm_dm_subdomains <bool> - use subdomains defined by the `DM` with `DMCreateDomainDecomposition()`
1724d709ea83SDmitry Karpeev 
1725d709ea83SDmitry Karpeev   Level: intermediate
1726d709ea83SDmitry Karpeev 
1727f1580f4eSBarry Smith   Note:
1728f1580f4eSBarry Smith   `PCASMSetTotalSubdomains()` and `PCASMSetOverlap()` take precedence over `PCASMSetDMSubdomains()`,
1729d709ea83SDmitry Karpeev   so setting either of the first two effectively turns the latter off.
1730d709ea83SDmitry Karpeev 
173173ff1848SBarry Smith   Developer Note:
173273ff1848SBarry Smith   This should be `PCASMSetUseDMSubdomains()`, similarly for the options database key
173373ff1848SBarry Smith 
1734562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMGetDMSubdomains()`, `PCASMSetTotalSubdomains()`, `PCASMSetOverlap()`
1735db781477SPatrick Sanan           `PCASMCreateSubdomains2D()`, `PCASMSetLocalSubdomains()`, `PCASMGetLocalSubdomains()`
1736d709ea83SDmitry Karpeev @*/
PCASMSetDMSubdomains(PC pc,PetscBool flg)1737d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMSetDMSubdomains(PC pc, PetscBool flg)
1738d71ae5a4SJacob Faibussowitsch {
1739d709ea83SDmitry Karpeev   PC_ASM   *osm = (PC_ASM *)pc->data;
1740d709ea83SDmitry Karpeev   PetscBool match;
1741d709ea83SDmitry Karpeev 
1742d709ea83SDmitry Karpeev   PetscFunctionBegin;
1743d709ea83SDmitry Karpeev   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1744d709ea83SDmitry Karpeev   PetscValidLogicalCollectiveBool(pc, flg, 2);
174528b400f6SJacob Faibussowitsch   PetscCheck(!pc->setupcalled, ((PetscObject)pc)->comm, PETSC_ERR_ARG_WRONGSTATE, "Not for a setup PC.");
17469566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCASM, &match));
1747ad540459SPierre Jolivet   if (match) osm->dm_subdomains = flg;
17483ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1749d709ea83SDmitry Karpeev }
1750d709ea83SDmitry Karpeev 
1751d709ea83SDmitry Karpeev /*@
1752f1580f4eSBarry Smith   PCASMGetDMSubdomains - Returns flag indicating whether to use `DMCreateDomainDecomposition()` to define the subdomains, whenever possible.
1753f1580f4eSBarry Smith 
1754d709ea83SDmitry Karpeev   Not Collective
1755d709ea83SDmitry Karpeev 
1756d709ea83SDmitry Karpeev   Input Parameter:
1757d709ea83SDmitry Karpeev . pc - the preconditioner
1758d709ea83SDmitry Karpeev 
1759d709ea83SDmitry Karpeev   Output Parameter:
176020f4b53cSBarry Smith . flg - boolean indicating whether to use subdomains defined by the `DM`
1761d709ea83SDmitry Karpeev 
1762d709ea83SDmitry Karpeev   Level: intermediate
1763d709ea83SDmitry Karpeev 
176473ff1848SBarry Smith   Developer Note:
176573ff1848SBarry Smith   This should be `PCASMSetUseDMSubdomains()`
176673ff1848SBarry Smith 
1767562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMSetDMSubdomains()`, `PCASMSetTotalSubdomains()`, `PCASMSetOverlap()`
1768db781477SPatrick Sanan           `PCASMCreateSubdomains2D()`, `PCASMSetLocalSubdomains()`, `PCASMGetLocalSubdomains()`
1769d709ea83SDmitry Karpeev @*/
PCASMGetDMSubdomains(PC pc,PetscBool * flg)1770d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMGetDMSubdomains(PC pc, PetscBool *flg)
1771d71ae5a4SJacob Faibussowitsch {
1772d709ea83SDmitry Karpeev   PC_ASM   *osm = (PC_ASM *)pc->data;
1773d709ea83SDmitry Karpeev   PetscBool match;
1774d709ea83SDmitry Karpeev 
1775d709ea83SDmitry Karpeev   PetscFunctionBegin;
1776d709ea83SDmitry Karpeev   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
17774f572ea9SToby Isaac   PetscAssertPointer(flg, 2);
17789566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCASM, &match));
177956ea09ceSStefano Zampini   if (match) *flg = osm->dm_subdomains;
178056ea09ceSStefano Zampini   else *flg = PETSC_FALSE;
17813ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1782d709ea83SDmitry Karpeev }
178380ec0b7dSPatrick Sanan 
17845d83a8b1SBarry Smith /*@
1785f1580f4eSBarry Smith   PCASMGetSubMatType - Gets the matrix type used for `PCASM` subsolves, as a string.
178680ec0b7dSPatrick Sanan 
178780ec0b7dSPatrick Sanan   Not Collective
178880ec0b7dSPatrick Sanan 
178980ec0b7dSPatrick Sanan   Input Parameter:
1790f1580f4eSBarry Smith . pc - the `PC`
179180ec0b7dSPatrick Sanan 
179280ec0b7dSPatrick Sanan   Output Parameter:
1793feefa0e1SJacob Faibussowitsch . sub_mat_type - name of matrix type
179480ec0b7dSPatrick Sanan 
179580ec0b7dSPatrick Sanan   Level: advanced
179680ec0b7dSPatrick Sanan 
1797562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMSetSubMatType()`, `PCSetType()`, `VecSetType()`, `MatType`, `Mat`
179880ec0b7dSPatrick Sanan @*/
PCASMGetSubMatType(PC pc,MatType * sub_mat_type)1799d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMGetSubMatType(PC pc, MatType *sub_mat_type)
1800d71ae5a4SJacob Faibussowitsch {
180156ea09ceSStefano Zampini   PetscFunctionBegin;
180256ea09ceSStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1803cac4c232SBarry Smith   PetscTryMethod(pc, "PCASMGetSubMatType_C", (PC, MatType *), (pc, sub_mat_type));
18043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
180580ec0b7dSPatrick Sanan }
180680ec0b7dSPatrick Sanan 
18075d83a8b1SBarry Smith /*@
1808f1580f4eSBarry Smith   PCASMSetSubMatType - Set the type of matrix used for `PCASM` subsolves
180980ec0b7dSPatrick Sanan 
1810c3339decSBarry Smith   Collective
181180ec0b7dSPatrick Sanan 
181280ec0b7dSPatrick Sanan   Input Parameters:
1813f1580f4eSBarry Smith + pc           - the `PC` object
1814f1580f4eSBarry Smith - sub_mat_type - the `MatType`
181580ec0b7dSPatrick Sanan 
181680ec0b7dSPatrick Sanan   Options Database Key:
1817f1580f4eSBarry Smith . -pc_asm_sub_mat_type  <sub_mat_type> - Sets the matrix type used for subsolves, for example, seqaijviennacl.
1818f1580f4eSBarry Smith    If you specify a base name like aijviennacl, the corresponding sequential type is assumed.
181980ec0b7dSPatrick Sanan 
1820f1580f4eSBarry Smith   Note:
18212fe279fdSBarry Smith   See `MatType` for available types
182280ec0b7dSPatrick Sanan 
182380ec0b7dSPatrick Sanan   Level: advanced
182480ec0b7dSPatrick Sanan 
1825562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMGetSubMatType()`, `PCSetType()`, `VecSetType()`, `MatType`, `Mat`
182680ec0b7dSPatrick Sanan @*/
PCASMSetSubMatType(PC pc,MatType sub_mat_type)1827d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMSetSubMatType(PC pc, MatType sub_mat_type)
1828d71ae5a4SJacob Faibussowitsch {
182956ea09ceSStefano Zampini   PetscFunctionBegin;
183056ea09ceSStefano Zampini   PetscValidHeaderSpecific(pc, PC_CLASSID, 1);
1831cac4c232SBarry Smith   PetscTryMethod(pc, "PCASMSetSubMatType_C", (PC, MatType), (pc, sub_mat_type));
18323ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
183380ec0b7dSPatrick Sanan }
1834