14b9ad928SBarry Smith /* 24b9ad928SBarry Smith This file defines an additive Schwarz preconditioner for any Mat implementation. 34b9ad928SBarry Smith 44b9ad928SBarry Smith Note that each processor may have any number of subdomains. But in order to 54b9ad928SBarry Smith deal easily with the VecScatter(), we treat each processor as if it has the 64b9ad928SBarry Smith same number of subdomains. 74b9ad928SBarry Smith 84b9ad928SBarry Smith n - total number of true subdomains on all processors 94b9ad928SBarry Smith n_local_true - actual number of subdomains on this processor 104b9ad928SBarry Smith n_local = maximum over all processors of n_local_true 114b9ad928SBarry Smith */ 124b9ad928SBarry Smith 1334eca32bSPierre Jolivet #include <petsc/private/pcasmimpl.h> /*I "petscpc.h" I*/ 1446233b44SBarry Smith #include <petsc/private/matimpl.h> 154b9ad928SBarry Smith 16d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCView_ASM(PC pc, PetscViewer viewer) 17d71ae5a4SJacob Faibussowitsch { 1892bb6962SLisandro Dalcin PC_ASM *osm = (PC_ASM *)pc->data; 1913f74950SBarry Smith PetscMPIInt rank; 204d219a6aSLisandro Dalcin PetscInt i, bsz; 21ace3abfcSBarry Smith PetscBool iascii, isstring; 224b9ad928SBarry Smith PetscViewer sviewer; 23ed155784SPierre Jolivet PetscViewerFormat format; 24ed155784SPierre Jolivet const char *prefix; 254b9ad928SBarry Smith 264b9ad928SBarry Smith PetscFunctionBegin; 279566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 289566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSTRING, &isstring)); 2932077d6dSBarry Smith if (iascii) { 303d03bb48SJed Brown char overlaps[256] = "user-defined overlap", blocks[256] = "total subdomain blocks not yet set"; 3163a3b9bcSJacob Faibussowitsch if (osm->overlap >= 0) PetscCall(PetscSNPrintf(overlaps, sizeof(overlaps), "amount of overlap = %" PetscInt_FMT, osm->overlap)); 3263a3b9bcSJacob Faibussowitsch if (osm->n > 0) PetscCall(PetscSNPrintf(blocks, sizeof(blocks), "total subdomain blocks = %" PetscInt_FMT, osm->n)); 339566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " %s, %s\n", blocks, overlaps)); 349566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " restriction/interpolation type - %s\n", PCASMTypes[osm->type])); 359566063dSJacob Faibussowitsch if (osm->dm_subdomains) PetscCall(PetscViewerASCIIPrintf(viewer, " Additive Schwarz: using DM to define subdomains\n")); 369566063dSJacob Faibussowitsch if (osm->loctype != PC_COMPOSITE_ADDITIVE) PetscCall(PetscViewerASCIIPrintf(viewer, " Additive Schwarz: local solve composition type - %s\n", PCCompositeTypes[osm->loctype])); 379566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank)); 389566063dSJacob Faibussowitsch PetscCall(PetscViewerGetFormat(viewer, &format)); 39ed155784SPierre Jolivet if (format != PETSC_VIEWER_ASCII_INFO_DETAIL) { 404d219a6aSLisandro Dalcin if (osm->ksp) { 419566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Local solver information for first block is in the following KSP and PC objects on rank 0:\n")); 429566063dSJacob Faibussowitsch PetscCall(PCGetOptionsPrefix(pc, &prefix)); 439566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Use -%sksp_view ::ascii_info_detail to display information for all blocks\n", prefix ? prefix : "")); 449566063dSJacob Faibussowitsch PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer)); 45dd400576SPatrick Sanan if (rank == 0) { 46b4025f61SBarry Smith PetscCall(PetscViewerASCIIPushTab(sviewer)); 479566063dSJacob Faibussowitsch PetscCall(KSPView(osm->ksp[0], sviewer)); 48b4025f61SBarry Smith PetscCall(PetscViewerASCIIPopTab(sviewer)); 494b9ad928SBarry Smith } 509566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer)); 514d219a6aSLisandro Dalcin } 524b9ad928SBarry Smith } else { 539566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 54835f2295SStefano Zampini PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " [%d] number of local blocks = %" PetscInt_FMT "\n", rank, osm->n_local_true)); 559566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 569566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Local solver information for each block is in the following KSP and PC objects:\n")); 579566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 589566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "- - - - - - - - - - - - - - - - - -\n")); 599566063dSJacob Faibussowitsch PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer)); 60dd2fa690SBarry Smith for (i = 0; i < osm->n_local_true; i++) { 619566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(osm->is[i], &bsz)); 62835f2295SStefano Zampini PetscCall(PetscViewerASCIIPrintf(sviewer, "[%d] local block number %" PetscInt_FMT ", size = %" PetscInt_FMT "\n", rank, i, bsz)); 639566063dSJacob Faibussowitsch PetscCall(KSPView(osm->ksp[i], sviewer)); 64b4025f61SBarry Smith PetscCall(PetscViewerASCIIPrintf(sviewer, "- - - - - - - - - - - - - - - - - -\n")); 654b9ad928SBarry Smith } 669566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer)); 679566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 689566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 694b9ad928SBarry Smith } 704b9ad928SBarry Smith } else if (isstring) { 7163a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerStringSPrintf(viewer, " blocks=%" PetscInt_FMT ", overlap=%" PetscInt_FMT ", type=%s", osm->n, osm->overlap, PCASMTypes[osm->type])); 729566063dSJacob Faibussowitsch PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer)); 739566063dSJacob Faibussowitsch if (osm->ksp) PetscCall(KSPView(osm->ksp[0], sviewer)); 749566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer)); 754b9ad928SBarry Smith } 763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 774b9ad928SBarry Smith } 784b9ad928SBarry Smith 79d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCASMPrintSubdomains(PC pc) 80d71ae5a4SJacob Faibussowitsch { 8192bb6962SLisandro Dalcin PC_ASM *osm = (PC_ASM *)pc->data; 8292bb6962SLisandro Dalcin const char *prefix; 8392bb6962SLisandro Dalcin char fname[PETSC_MAX_PATH_LEN + 1]; 84643535aeSDmitry Karpeev PetscViewer viewer, sviewer; 85643535aeSDmitry Karpeev char *s; 8692bb6962SLisandro Dalcin PetscInt i, j, nidx; 8792bb6962SLisandro Dalcin const PetscInt *idx; 88643535aeSDmitry Karpeev PetscMPIInt rank, size; 89846783a0SBarry Smith 9092bb6962SLisandro Dalcin PetscFunctionBegin; 919566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size)); 929566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank)); 939566063dSJacob Faibussowitsch PetscCall(PCGetOptionsPrefix(pc, &prefix)); 949566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, prefix, "-pc_asm_print_subdomains", fname, sizeof(fname), NULL)); 95c6a7a370SJeremy L Thompson if (fname[0] == 0) PetscCall(PetscStrncpy(fname, "stdout", sizeof(fname))); 969566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIOpen(PetscObjectComm((PetscObject)pc), fname, &viewer)); 97643535aeSDmitry Karpeev for (i = 0; i < osm->n_local; i++) { 98643535aeSDmitry Karpeev if (i < osm->n_local_true) { 999566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(osm->is[i], &nidx)); 1009566063dSJacob Faibussowitsch PetscCall(ISGetIndices(osm->is[i], &idx)); 101643535aeSDmitry Karpeev /* Print to a string viewer; no more than 15 characters per index plus 512 char for the header.*/ 10236a9e3b9SBarry Smith #define len 16 * (nidx + 1) + 512 1039566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(len, &s)); 1049566063dSJacob Faibussowitsch PetscCall(PetscViewerStringOpen(PETSC_COMM_SELF, s, len, &sviewer)); 10536a9e3b9SBarry Smith #undef len 10663a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerStringSPrintf(sviewer, "[%d:%d] Subdomain %" PetscInt_FMT " with overlap:\n", rank, size, i)); 10748a46eb9SPierre Jolivet for (j = 0; j < nidx; j++) PetscCall(PetscViewerStringSPrintf(sviewer, "%" PetscInt_FMT " ", idx[j])); 1089566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(osm->is[i], &idx)); 1099566063dSJacob Faibussowitsch PetscCall(PetscViewerStringSPrintf(sviewer, "\n")); 1109566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&sviewer)); 1119566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 11263a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "%s", s)); 1139566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 1149566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 1159566063dSJacob Faibussowitsch PetscCall(PetscFree(s)); 1162b691e39SMatthew Knepley if (osm->is_local) { 117643535aeSDmitry Karpeev /* Print to a string viewer; no more than 15 characters per index plus 512 char for the header.*/ 11836a9e3b9SBarry Smith #define len 16 * (nidx + 1) + 512 1199566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(len, &s)); 1209566063dSJacob Faibussowitsch PetscCall(PetscViewerStringOpen(PETSC_COMM_SELF, s, len, &sviewer)); 12136a9e3b9SBarry Smith #undef len 12263a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerStringSPrintf(sviewer, "[%d:%d] Subdomain %" PetscInt_FMT " without overlap:\n", rank, size, i)); 1239566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(osm->is_local[i], &nidx)); 1249566063dSJacob Faibussowitsch PetscCall(ISGetIndices(osm->is_local[i], &idx)); 12548a46eb9SPierre Jolivet for (j = 0; j < nidx; j++) PetscCall(PetscViewerStringSPrintf(sviewer, "%" PetscInt_FMT " ", idx[j])); 1269566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(osm->is_local[i], &idx)); 1279566063dSJacob Faibussowitsch PetscCall(PetscViewerStringSPrintf(sviewer, "\n")); 1289566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&sviewer)); 1299566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 13063a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "%s", s)); 1319566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 1329566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 1339566063dSJacob Faibussowitsch PetscCall(PetscFree(s)); 134643535aeSDmitry Karpeev } 1352fa5cd67SKarl Rupp } else { 136643535aeSDmitry Karpeev /* Participate in collective viewer calls. */ 1379566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 1389566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 1399566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 140643535aeSDmitry Karpeev /* Assume either all ranks have is_local or none do. */ 141643535aeSDmitry Karpeev if (osm->is_local) { 1429566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 1439566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 1449566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 145643535aeSDmitry Karpeev } 146fdc48646SDmitry Karpeev } 14792bb6962SLisandro Dalcin } 1489566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 1499566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 1503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 15192bb6962SLisandro Dalcin } 15292bb6962SLisandro Dalcin 153d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_ASM(PC pc) 154d71ae5a4SJacob Faibussowitsch { 1554b9ad928SBarry Smith PC_ASM *osm = (PC_ASM *)pc->data; 15657501b6eSBarry Smith PetscBool flg; 15787e86712SBarry Smith PetscInt i, m, m_local; 1584b9ad928SBarry Smith MatReuse scall = MAT_REUSE_MATRIX; 1594b9ad928SBarry Smith IS isl; 1604b9ad928SBarry Smith KSP ksp; 1614b9ad928SBarry Smith PC subpc; 1622dcb1b2aSMatthew Knepley const char *prefix, *pprefix; 16323ce1328SBarry Smith Vec vec; 1640298fd71SBarry Smith DM *domain_dm = NULL; 1654849c82aSBarry Smith MatNullSpace *nullsp = NULL; 1664b9ad928SBarry Smith 1674b9ad928SBarry Smith PetscFunctionBegin; 1684b9ad928SBarry Smith if (!pc->setupcalled) { 169265a2120SBarry Smith PetscInt m; 17092bb6962SLisandro Dalcin 1712b837212SDmitry Karpeev /* Note: if subdomains have been set either via PCASMSetTotalSubdomains() or via PCASMSetLocalSubdomains(), osm->n_local_true will not be PETSC_DECIDE */ 1722b837212SDmitry Karpeev if (osm->n_local_true == PETSC_DECIDE) { 17369ca1f37SDmitry Karpeev /* no subdomains given */ 17465db9045SDmitry Karpeev /* try pc->dm first, if allowed */ 175d709ea83SDmitry Karpeev if (osm->dm_subdomains && pc->dm) { 176feb221f8SDmitry Karpeev PetscInt num_domains, d; 177feb221f8SDmitry Karpeev char **domain_names; 1788d4ac253SDmitry Karpeev IS *inner_domain_is, *outer_domain_is; 1799566063dSJacob Faibussowitsch PetscCall(DMCreateDomainDecomposition(pc->dm, &num_domains, &domain_names, &inner_domain_is, &outer_domain_is, &domain_dm)); 180704f0395SPatrick Sanan osm->overlap = -1; /* We do not want to increase the overlap of the IS. 181704f0395SPatrick Sanan A future improvement of this code might allow one to use 182704f0395SPatrick Sanan DM-defined subdomains and also increase the overlap, 183704f0395SPatrick Sanan but that is not currently supported */ 1841baa6e33SBarry Smith if (num_domains) PetscCall(PCASMSetLocalSubdomains(pc, num_domains, outer_domain_is, inner_domain_is)); 185feb221f8SDmitry Karpeev for (d = 0; d < num_domains; ++d) { 1869566063dSJacob Faibussowitsch if (domain_names) PetscCall(PetscFree(domain_names[d])); 1879566063dSJacob Faibussowitsch if (inner_domain_is) PetscCall(ISDestroy(&inner_domain_is[d])); 1889566063dSJacob Faibussowitsch if (outer_domain_is) PetscCall(ISDestroy(&outer_domain_is[d])); 189feb221f8SDmitry Karpeev } 1909566063dSJacob Faibussowitsch PetscCall(PetscFree(domain_names)); 1919566063dSJacob Faibussowitsch PetscCall(PetscFree(inner_domain_is)); 1929566063dSJacob Faibussowitsch PetscCall(PetscFree(outer_domain_is)); 193feb221f8SDmitry Karpeev } 1942b837212SDmitry Karpeev if (osm->n_local_true == PETSC_DECIDE) { 19569ca1f37SDmitry Karpeev /* still no subdomains; use one subdomain per processor */ 1962b837212SDmitry Karpeev osm->n_local_true = 1; 197feb221f8SDmitry Karpeev } 1982b837212SDmitry Karpeev } 1992b837212SDmitry Karpeev { /* determine the global and max number of subdomains */ 2009371c9d4SSatish Balay struct { 2019371c9d4SSatish Balay PetscInt max, sum; 2029371c9d4SSatish Balay } inwork, outwork; 203c10200c1SHong Zhang PetscMPIInt size; 204c10200c1SHong Zhang 2056ac3741eSJed Brown inwork.max = osm->n_local_true; 2066ac3741eSJed Brown inwork.sum = osm->n_local_true; 207462c564dSBarry Smith PetscCallMPI(MPIU_Allreduce(&inwork, &outwork, 1, MPIU_2INT, MPIU_MAXSUM_OP, PetscObjectComm((PetscObject)pc))); 2086ac3741eSJed Brown osm->n_local = outwork.max; 2096ac3741eSJed Brown osm->n = outwork.sum; 210c10200c1SHong Zhang 2119566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size)); 212c10200c1SHong Zhang if (outwork.max == 1 && outwork.sum == size) { 2137dae84e0SHong Zhang /* osm->n_local_true = 1 on all processes, set this option may enable use of optimized MatCreateSubMatrices() implementation */ 2149566063dSJacob Faibussowitsch PetscCall(MatSetOption(pc->pmat, MAT_SUBMAT_SINGLEIS, PETSC_TRUE)); 215c10200c1SHong Zhang } 2164b9ad928SBarry Smith } 21778904715SLisandro Dalcin if (!osm->is) { /* create the index sets */ 2189566063dSJacob Faibussowitsch PetscCall(PCASMCreateSubdomains(pc->pmat, osm->n_local_true, &osm->is)); 2194b9ad928SBarry Smith } 220f5234e35SJed Brown if (osm->n_local_true > 1 && !osm->is_local) { 2219566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(osm->n_local_true, &osm->is_local)); 222f5234e35SJed Brown for (i = 0; i < osm->n_local_true; i++) { 223f5234e35SJed Brown if (osm->overlap > 0) { /* With positive overlap, osm->is[i] will be modified */ 2249566063dSJacob Faibussowitsch PetscCall(ISDuplicate(osm->is[i], &osm->is_local[i])); 2259566063dSJacob Faibussowitsch PetscCall(ISCopy(osm->is[i], osm->is_local[i])); 226f5234e35SJed Brown } else { 2279566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)osm->is[i])); 228f5234e35SJed Brown osm->is_local[i] = osm->is[i]; 229f5234e35SJed Brown } 230f5234e35SJed Brown } 231f5234e35SJed Brown } 2329566063dSJacob Faibussowitsch PetscCall(PCGetOptionsPrefix(pc, &prefix)); 2333d03bb48SJed Brown if (osm->overlap > 0) { 2344b9ad928SBarry Smith /* Extend the "overlapping" regions by a number of steps */ 2359566063dSJacob Faibussowitsch PetscCall(MatIncreaseOverlap(pc->pmat, osm->n_local_true, osm->is, osm->overlap)); 2363d03bb48SJed Brown } 2376ed231c7SMatthew Knepley if (osm->sort_indices) { 23892bb6962SLisandro Dalcin for (i = 0; i < osm->n_local_true; i++) { 2399566063dSJacob Faibussowitsch PetscCall(ISSort(osm->is[i])); 24048a46eb9SPierre Jolivet if (osm->is_local) PetscCall(ISSort(osm->is_local[i])); 2414b9ad928SBarry Smith } 2426ed231c7SMatthew Knepley } 24398d3e228SPierre Jolivet flg = PETSC_FALSE; 24498d3e228SPierre Jolivet PetscCall(PetscOptionsHasName(NULL, prefix, "-pc_asm_print_subdomains", &flg)); 24598d3e228SPierre Jolivet if (flg) PetscCall(PCASMPrintSubdomains(pc)); 246f6991133SBarry Smith if (!osm->ksp) { 24778904715SLisandro Dalcin /* Create the local solvers */ 2489566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(osm->n_local_true, &osm->ksp)); 24948a46eb9SPierre Jolivet if (domain_dm) PetscCall(PetscInfo(pc, "Setting up ASM subproblems using the embedded DM\n")); 25092bb6962SLisandro Dalcin for (i = 0; i < osm->n_local_true; i++) { 2519566063dSJacob Faibussowitsch PetscCall(KSPCreate(PETSC_COMM_SELF, &ksp)); 2523821be0aSBarry Smith PetscCall(KSPSetNestLevel(ksp, pc->kspnestlevel)); 2539566063dSJacob Faibussowitsch PetscCall(KSPSetErrorIfNotConverged(ksp, pc->erroriffailure)); 2549566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)ksp, (PetscObject)pc, 1)); 2559566063dSJacob Faibussowitsch PetscCall(KSPSetType(ksp, KSPPREONLY)); 2569566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp, &subpc)); 2579566063dSJacob Faibussowitsch PetscCall(PCGetOptionsPrefix(pc, &prefix)); 2589566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(ksp, prefix)); 2599566063dSJacob Faibussowitsch PetscCall(KSPAppendOptionsPrefix(ksp, "sub_")); 260feb221f8SDmitry Karpeev if (domain_dm) { 2619566063dSJacob Faibussowitsch PetscCall(KSPSetDM(ksp, domain_dm[i])); 2629566063dSJacob Faibussowitsch PetscCall(KSPSetDMActive(ksp, PETSC_FALSE)); 2639566063dSJacob Faibussowitsch PetscCall(DMDestroy(&domain_dm[i])); 264feb221f8SDmitry Karpeev } 2654b9ad928SBarry Smith osm->ksp[i] = ksp; 2664b9ad928SBarry Smith } 2671baa6e33SBarry Smith if (domain_dm) PetscCall(PetscFree(domain_dm)); 268f6991133SBarry Smith } 2691dd8081eSeaulisa 2709566063dSJacob Faibussowitsch PetscCall(ISConcatenate(PETSC_COMM_SELF, osm->n_local_true, osm->is, &osm->lis)); 2719566063dSJacob Faibussowitsch PetscCall(ISSortRemoveDups(osm->lis)); 2729566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(osm->lis, &m)); 2731dd8081eSeaulisa 2744b9ad928SBarry Smith scall = MAT_INITIAL_MATRIX; 2754b9ad928SBarry Smith } else { 2764b9ad928SBarry Smith /* 2774b9ad928SBarry Smith Destroy the blocks from the previous iteration 2784b9ad928SBarry Smith */ 279e09e08ccSBarry Smith if (pc->flag == DIFFERENT_NONZERO_PATTERN) { 2804849c82aSBarry Smith PetscCall(MatGetNullSpaces(osm->n_local_true, osm->pmat, &nullsp)); 2819566063dSJacob Faibussowitsch PetscCall(MatDestroyMatrices(osm->n_local_true, &osm->pmat)); 2824b9ad928SBarry Smith scall = MAT_INITIAL_MATRIX; 2834b9ad928SBarry Smith } 2844b9ad928SBarry Smith } 2854b9ad928SBarry Smith 286b58cb649SBarry Smith /* Destroy previous submatrices of a different type than pc->pmat since MAT_REUSE_MATRIX won't work in that case */ 28782b5ce2aSStefano Zampini if (scall == MAT_REUSE_MATRIX && osm->sub_mat_type) { 2884849c82aSBarry Smith PetscCall(MatGetNullSpaces(osm->n_local_true, osm->pmat, &nullsp)); 28948a46eb9SPierre Jolivet if (osm->n_local_true > 0) PetscCall(MatDestroySubMatrices(osm->n_local_true, &osm->pmat)); 290b58cb649SBarry Smith scall = MAT_INITIAL_MATRIX; 291b58cb649SBarry Smith } 292b58cb649SBarry Smith 29392bb6962SLisandro Dalcin /* 29492bb6962SLisandro Dalcin Extract out the submatrices 29592bb6962SLisandro Dalcin */ 2969566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(pc->pmat, osm->n_local_true, osm->is, osm->is, scall, &osm->pmat)); 29792bb6962SLisandro Dalcin if (scall == MAT_INITIAL_MATRIX) { 2989566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptionsPrefix((PetscObject)pc->pmat, &pprefix)); 299aa624791SPierre Jolivet for (i = 0; i < osm->n_local_true; i++) PetscCall(PetscObjectSetOptionsPrefix((PetscObject)osm->pmat[i], pprefix)); 3004849c82aSBarry Smith if (nullsp) PetscCall(MatRestoreNullSpaces(osm->n_local_true, osm->pmat, &nullsp)); 30192bb6962SLisandro Dalcin } 30280ec0b7dSPatrick Sanan 30380ec0b7dSPatrick Sanan /* Convert the types of the submatrices (if needbe) */ 30480ec0b7dSPatrick Sanan if (osm->sub_mat_type) { 305f4f49eeaSPierre Jolivet for (i = 0; i < osm->n_local_true; i++) PetscCall(MatConvert(osm->pmat[i], osm->sub_mat_type, MAT_INPLACE_MATRIX, &osm->pmat[i])); 30680ec0b7dSPatrick Sanan } 30780ec0b7dSPatrick Sanan 30880ec0b7dSPatrick Sanan if (!pc->setupcalled) { 30956ea09ceSStefano Zampini VecType vtype; 31056ea09ceSStefano Zampini 31180ec0b7dSPatrick Sanan /* Create the local work vectors (from the local matrices) and scatter contexts */ 3129566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(pc->pmat, &vec, NULL)); 3135b3c0d42Seaulisa 314677c7726SPierre Jolivet PetscCheck(!osm->is_local || osm->n_local_true == 1 || (osm->type != PC_ASM_INTERPOLATE && osm->type != PC_ASM_NONE), PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Cannot use interpolate or none PCASMType if is_local was provided to PCASMSetLocalSubdomains() with more than a single subdomain"); 315677c7726SPierre Jolivet if (osm->is_local && osm->type != PC_ASM_BASIC && osm->loctype == PC_COMPOSITE_ADDITIVE) PetscCall(PetscMalloc1(osm->n_local_true, &osm->lprolongation)); 3169566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(osm->n_local_true, &osm->lrestriction)); 3179566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(osm->n_local_true, &osm->x)); 3189566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(osm->n_local_true, &osm->y)); 319347574c9Seaulisa 3209566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(osm->lis, &m)); 3219566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, m, 0, 1, &isl)); 3229566063dSJacob Faibussowitsch PetscCall(MatGetVecType(osm->pmat[0], &vtype)); 3239566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_SELF, &osm->lx)); 3249566063dSJacob Faibussowitsch PetscCall(VecSetSizes(osm->lx, m, m)); 3259566063dSJacob Faibussowitsch PetscCall(VecSetType(osm->lx, vtype)); 3269566063dSJacob Faibussowitsch PetscCall(VecDuplicate(osm->lx, &osm->ly)); 3279566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(vec, osm->lis, osm->lx, isl, &osm->restriction)); 3289566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isl)); 329347574c9Seaulisa 33080ec0b7dSPatrick Sanan for (i = 0; i < osm->n_local_true; ++i) { 3315b3c0d42Seaulisa ISLocalToGlobalMapping ltog; 3325b3c0d42Seaulisa IS isll; 3335b3c0d42Seaulisa const PetscInt *idx_is; 3345b3c0d42Seaulisa PetscInt *idx_lis, nout; 3355b3c0d42Seaulisa 3369566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(osm->is[i], &m)); 3379566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(osm->pmat[i], &osm->x[i], NULL)); 3389566063dSJacob Faibussowitsch PetscCall(VecDuplicate(osm->x[i], &osm->y[i])); 33955817e79SBarry Smith 340b0de9f23SBarry Smith /* generate a scatter from ly to y[i] picking all the overlapping is[i] entries */ 3419566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreateIS(osm->lis, <og)); 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(<og)); 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], <og)); 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(<og)); 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, <og)); 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(<og)); 3747827d75bSBarry Smith PetscCheck(nout == m_local, PETSC_COMM_SELF, PETSC_ERR_PLIB, "is_local not a subset of lis"); 3759566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, m_local, idx2, PETSC_OWN_POINTER, &isll_local)); 376d8b3b5e3Seaulisa 3779566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(osm->is_local[i], &idx_local)); 3789566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(osm->y[i], isll, osm->ly, isll_local, &osm->lprolongation[i])); 379d8b3b5e3Seaulisa 3809566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isll)); 3819566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isll_local)); 382d8b3b5e3Seaulisa } 38380ec0b7dSPatrick Sanan } 3849566063dSJacob Faibussowitsch PetscCall(VecDestroy(&vec)); 38580ec0b7dSPatrick Sanan } 38680ec0b7dSPatrick Sanan 387fb745f2cSMatthew G. Knepley if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE) { 388235cc4ceSMatthew G. Knepley IS *cis; 389235cc4ceSMatthew G. Knepley PetscInt c; 390235cc4ceSMatthew G. Knepley 3919566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(osm->n_local_true, &cis)); 392235cc4ceSMatthew G. Knepley for (c = 0; c < osm->n_local_true; ++c) cis[c] = osm->lis; 3939566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(pc->pmat, osm->n_local_true, osm->is, cis, scall, &osm->lmats)); 3949566063dSJacob Faibussowitsch PetscCall(PetscFree(cis)); 395fb745f2cSMatthew G. Knepley } 3964b9ad928SBarry Smith 3974b9ad928SBarry Smith /* Return control to the user so that the submatrices can be modified (e.g., to apply 3984b9ad928SBarry Smith different boundary conditions for the submatrices than for the global problem) */ 3999566063dSJacob Faibussowitsch PetscCall(PCModifySubMatrices(pc, osm->n_local_true, osm->is, osm->is, osm->pmat, pc->modifysubmatricesP)); 4004b9ad928SBarry Smith 40192bb6962SLisandro Dalcin /* 40292bb6962SLisandro Dalcin Loop over subdomains putting them into local ksp 40392bb6962SLisandro Dalcin */ 4049566063dSJacob Faibussowitsch PetscCall(KSPGetOptionsPrefix(osm->ksp[0], &prefix)); 40592bb6962SLisandro Dalcin for (i = 0; i < osm->n_local_true; i++) { 4069566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(osm->ksp[i], osm->pmat[i], osm->pmat[i])); 4079566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(osm->pmat[i], prefix)); 40848a46eb9SPierre Jolivet if (!pc->setupcalled) PetscCall(KSPSetFromOptions(osm->ksp[i])); 40992bb6962SLisandro Dalcin } 4103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4114b9ad928SBarry Smith } 4124b9ad928SBarry Smith 413d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUpOnBlocks_ASM(PC pc) 414d71ae5a4SJacob Faibussowitsch { 4154b9ad928SBarry Smith PC_ASM *osm = (PC_ASM *)pc->data; 41613f74950SBarry Smith PetscInt i; 417539c167fSBarry Smith KSPConvergedReason reason; 4184b9ad928SBarry Smith 4194b9ad928SBarry Smith PetscFunctionBegin; 4204b9ad928SBarry Smith for (i = 0; i < osm->n_local_true; i++) { 4219566063dSJacob Faibussowitsch PetscCall(KSPSetUp(osm->ksp[i])); 4229566063dSJacob Faibussowitsch PetscCall(KSPGetConvergedReason(osm->ksp[i], &reason)); 423ad540459SPierre Jolivet if (reason == KSP_DIVERGED_PC_FAILED) pc->failedreason = PC_SUBPC_ERROR; 4244b9ad928SBarry Smith } 4253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4264b9ad928SBarry Smith } 4274b9ad928SBarry Smith 428d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_ASM(PC pc, Vec x, Vec y) 429d71ae5a4SJacob Faibussowitsch { 4304b9ad928SBarry Smith PC_ASM *osm = (PC_ASM *)pc->data; 4311dd8081eSeaulisa PetscInt i, n_local_true = osm->n_local_true; 4324b9ad928SBarry Smith ScatterMode forward = SCATTER_FORWARD, reverse = SCATTER_REVERSE; 4334b9ad928SBarry Smith 4344b9ad928SBarry Smith PetscFunctionBegin; 4354b9ad928SBarry Smith /* 43648e38a8aSPierre Jolivet support for limiting the restriction or interpolation to only local 4374b9ad928SBarry Smith subdomain values (leaving the other values 0). 4384b9ad928SBarry Smith */ 4394b9ad928SBarry Smith if (!(osm->type & PC_ASM_RESTRICT)) { 4404b9ad928SBarry Smith forward = SCATTER_FORWARD_LOCAL; 4414b9ad928SBarry Smith /* have to zero the work RHS since scatter may leave some slots empty */ 4429566063dSJacob Faibussowitsch PetscCall(VecSet(osm->lx, 0.0)); 4434b9ad928SBarry Smith } 444ad540459SPierre Jolivet if (!(osm->type & PC_ASM_INTERPOLATE)) reverse = SCATTER_REVERSE_LOCAL; 4454b9ad928SBarry Smith 4460fdf79fbSJacob Faibussowitsch PetscCheck(osm->loctype == PC_COMPOSITE_MULTIPLICATIVE || osm->loctype == PC_COMPOSITE_ADDITIVE, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Invalid local composition type: %s", PCCompositeTypes[osm->loctype]); 447b0de9f23SBarry Smith /* zero the global and the local solutions */ 4489566063dSJacob Faibussowitsch PetscCall(VecSet(y, 0.0)); 4499566063dSJacob Faibussowitsch PetscCall(VecSet(osm->ly, 0.0)); 450347574c9Seaulisa 45148e38a8aSPierre Jolivet /* copy the global RHS to local RHS including the ghost nodes */ 4529566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->restriction, x, osm->lx, INSERT_VALUES, forward)); 4539566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->restriction, x, osm->lx, INSERT_VALUES, forward)); 454347574c9Seaulisa 45548e38a8aSPierre Jolivet /* restrict local RHS to the overlapping 0-block RHS */ 4569566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward)); 4579566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward)); 458d8b3b5e3Seaulisa 45912cd4985SMatthew G. Knepley /* do the local solves */ 46012cd4985SMatthew G. Knepley for (i = 0; i < n_local_true; ++i) { 461b0de9f23SBarry Smith /* solve the overlapping i-block */ 4629566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PC_ApplyOnBlocks, osm->ksp[i], osm->x[i], osm->y[i], 0)); 4639566063dSJacob Faibussowitsch PetscCall(KSPSolve(osm->ksp[i], osm->x[i], osm->y[i])); 4649566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(osm->ksp[i], pc, osm->y[i])); 4659566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PC_ApplyOnBlocks, osm->ksp[i], osm->x[i], osm->y[i], 0)); 466d8b3b5e3Seaulisa 467677c7726SPierre Jolivet if (osm->lprolongation && osm->type != PC_ASM_INTERPOLATE) { /* interpolate the non-overlapping i-block solution to the local solution (only for restrictive additive) */ 4689566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward)); 4699566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward)); 47048e38a8aSPierre Jolivet } else { /* interpolate the overlapping i-block solution to the local solution */ 4719566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse)); 4729566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse)); 473d8b3b5e3Seaulisa } 474347574c9Seaulisa 475347574c9Seaulisa if (i < n_local_true - 1) { 47648e38a8aSPierre Jolivet /* restrict local RHS to the overlapping (i+1)-block RHS */ 4779566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->lrestriction[i + 1], osm->lx, osm->x[i + 1], INSERT_VALUES, forward)); 4789566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->lrestriction[i + 1], osm->lx, osm->x[i + 1], INSERT_VALUES, forward)); 479347574c9Seaulisa 480347574c9Seaulisa if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE) { 4818966356dSPierre Jolivet /* update the overlapping (i+1)-block RHS using the current local solution */ 4829566063dSJacob Faibussowitsch PetscCall(MatMult(osm->lmats[i + 1], osm->ly, osm->y[i + 1])); 4839566063dSJacob Faibussowitsch PetscCall(VecAXPBY(osm->x[i + 1], -1., 1., osm->y[i + 1])); 4847c3d802fSMatthew G. Knepley } 48512cd4985SMatthew G. Knepley } 48612cd4985SMatthew G. Knepley } 48748e38a8aSPierre Jolivet /* add the local solution to the global solution including the ghost nodes */ 4889566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->restriction, osm->ly, y, ADD_VALUES, reverse)); 4899566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->restriction, osm->ly, y, ADD_VALUES, reverse)); 4903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4914b9ad928SBarry Smith } 4924b9ad928SBarry Smith 493d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMatApply_ASM(PC pc, Mat X, Mat Y) 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 */ 50748e38a8aSPierre Jolivet if (!(osm->type & PC_ASM_RESTRICT)) { 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 } 512ad540459SPierre Jolivet if (!(osm->type & PC_ASM_INTERPOLATE)) 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 */ 5379566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PC_ApplyOnBlocks, osm->ksp[0], Z, W, 0)); 5389566063dSJacob Faibussowitsch PetscCall(KSPMatSolve(osm->ksp[0], Z, W)); 5399566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(osm->ksp[0], pc, NULL)); 5409566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PC_ApplyOnBlocks, osm->ksp[0], Z, W, 0)); 5419566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Z)); 54248e38a8aSPierre Jolivet 54348e38a8aSPierre Jolivet for (i = 0; i < N; ++i) { 5449566063dSJacob Faibussowitsch PetscCall(VecSet(osm->ly, 0.0)); 5459566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumnVecRead(W, i, &x)); 546677c7726SPierre Jolivet if (osm->lprolongation && osm->type != PC_ASM_INTERPOLATE) { /* interpolate the non-overlapping 0-block solution to the local solution (only for restrictive additive) */ 5479566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->lprolongation[0], x, osm->ly, ADD_VALUES, forward)); 5489566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->lprolongation[0], x, osm->ly, ADD_VALUES, forward)); 54948e38a8aSPierre Jolivet } else { /* interpolate the overlapping 0-block solution to the local solution */ 5509566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->lrestriction[0], x, osm->ly, ADD_VALUES, reverse)); 5519566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->lrestriction[0], x, osm->ly, ADD_VALUES, reverse)); 55248e38a8aSPierre Jolivet } 5539566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumnVecRead(W, i, &x)); 55448e38a8aSPierre Jolivet 5559566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumnVecWrite(Y, i, &x)); 55648e38a8aSPierre Jolivet /* add the local solution to the global solution including the ghost nodes */ 5579566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->restriction, osm->ly, x, ADD_VALUES, reverse)); 5589566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->restriction, osm->ly, x, ADD_VALUES, reverse)); 5599566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumnVecWrite(Y, i, &x)); 56048e38a8aSPierre Jolivet } 5619566063dSJacob Faibussowitsch PetscCall(MatDestroy(&W)); 5623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 56348e38a8aSPierre Jolivet } 56448e38a8aSPierre Jolivet 565d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplyTranspose_ASM(PC pc, Vec x, Vec y) 566d71ae5a4SJacob Faibussowitsch { 5674b9ad928SBarry Smith PC_ASM *osm = (PC_ASM *)pc->data; 5681dd8081eSeaulisa PetscInt i, n_local_true = osm->n_local_true; 5694b9ad928SBarry Smith ScatterMode forward = SCATTER_FORWARD, reverse = SCATTER_REVERSE; 5704b9ad928SBarry Smith 5714b9ad928SBarry Smith PetscFunctionBegin; 572*6040bc1fSPierre Jolivet PetscCheck(osm->n_local_true <= 1 || osm->loctype == PC_COMPOSITE_ADDITIVE, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Not yet implemented"); 5734b9ad928SBarry Smith /* 5744b9ad928SBarry Smith Support for limiting the restriction or interpolation to only local 5754b9ad928SBarry Smith subdomain values (leaving the other values 0). 5764b9ad928SBarry Smith 5774b9ad928SBarry Smith Note: these are reversed from the PCApply_ASM() because we are applying the 5784b9ad928SBarry Smith transpose of the three terms 5794b9ad928SBarry Smith */ 580d8b3b5e3Seaulisa 5814b9ad928SBarry Smith if (!(osm->type & PC_ASM_INTERPOLATE)) { 5824b9ad928SBarry Smith forward = SCATTER_FORWARD_LOCAL; 5834b9ad928SBarry Smith /* have to zero the work RHS since scatter may leave some slots empty */ 5849566063dSJacob Faibussowitsch PetscCall(VecSet(osm->lx, 0.0)); 5854b9ad928SBarry Smith } 5862fa5cd67SKarl Rupp if (!(osm->type & PC_ASM_RESTRICT)) reverse = SCATTER_REVERSE_LOCAL; 5874b9ad928SBarry Smith 588b0de9f23SBarry Smith /* zero the global and the local solutions */ 5899566063dSJacob Faibussowitsch PetscCall(VecSet(y, 0.0)); 5909566063dSJacob Faibussowitsch PetscCall(VecSet(osm->ly, 0.0)); 591d8b3b5e3Seaulisa 592b0de9f23SBarry Smith /* Copy the global RHS to local RHS including the ghost nodes */ 5939566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->restriction, x, osm->lx, INSERT_VALUES, forward)); 5949566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->restriction, x, osm->lx, INSERT_VALUES, forward)); 595d8b3b5e3Seaulisa 596b0de9f23SBarry Smith /* Restrict local RHS to the overlapping 0-block RHS */ 5979566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward)); 5989566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward)); 599d8b3b5e3Seaulisa 6004b9ad928SBarry Smith /* do the local solves */ 601d8b3b5e3Seaulisa for (i = 0; i < n_local_true; ++i) { 602b0de9f23SBarry Smith /* solve the overlapping i-block */ 603*6040bc1fSPierre Jolivet PetscCall(PetscLogEventBegin(PC_ApplyTransposeOnBlocks, osm->ksp[i], osm->x[i], osm->y[i], 0)); 6049566063dSJacob Faibussowitsch PetscCall(KSPSolveTranspose(osm->ksp[i], osm->x[i], osm->y[i])); 6059566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(osm->ksp[i], pc, osm->y[i])); 606*6040bc1fSPierre Jolivet PetscCall(PetscLogEventEnd(PC_ApplyTransposeOnBlocks, osm->ksp[i], osm->x[i], osm->y[i], 0)); 607d8b3b5e3Seaulisa 608a681a0f1SPierre Jolivet if (osm->lprolongation && osm->type != PC_ASM_RESTRICT) { /* interpolate the non-overlapping i-block solution to the local solution */ 6099566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward)); 6109566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward)); 611910cf402Sprj- } else { /* interpolate the overlapping i-block solution to the local solution */ 6129566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse)); 6139566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse)); 6144b9ad928SBarry Smith } 615d8b3b5e3Seaulisa 616d8b3b5e3Seaulisa if (i < n_local_true - 1) { 617b0de9f23SBarry Smith /* Restrict local RHS to the overlapping (i+1)-block RHS */ 6189566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->lrestriction[i + 1], osm->lx, osm->x[i + 1], INSERT_VALUES, forward)); 6199566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->lrestriction[i + 1], osm->lx, osm->x[i + 1], INSERT_VALUES, forward)); 6204b9ad928SBarry Smith } 6214b9ad928SBarry Smith } 622b0de9f23SBarry Smith /* Add the local solution to the global solution including the ghost nodes */ 6239566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->restriction, osm->ly, y, ADD_VALUES, reverse)); 6249566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->restriction, osm->ly, y, ADD_VALUES, reverse)); 6253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6264b9ad928SBarry Smith } 6274b9ad928SBarry Smith 628d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCReset_ASM(PC pc) 629d71ae5a4SJacob Faibussowitsch { 6304b9ad928SBarry Smith PC_ASM *osm = (PC_ASM *)pc->data; 63113f74950SBarry Smith PetscInt i; 6324b9ad928SBarry Smith 6334b9ad928SBarry Smith PetscFunctionBegin; 63492bb6962SLisandro Dalcin if (osm->ksp) { 63548a46eb9SPierre Jolivet for (i = 0; i < osm->n_local_true; i++) PetscCall(KSPReset(osm->ksp[i])); 63692bb6962SLisandro Dalcin } 637e09e08ccSBarry Smith if (osm->pmat) { 63848a46eb9SPierre Jolivet if (osm->n_local_true > 0) PetscCall(MatDestroySubMatrices(osm->n_local_true, &osm->pmat)); 63992bb6962SLisandro Dalcin } 6401dd8081eSeaulisa if (osm->lrestriction) { 6419566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&osm->restriction)); 6421dd8081eSeaulisa for (i = 0; i < osm->n_local_true; i++) { 6439566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&osm->lrestriction[i])); 6449566063dSJacob Faibussowitsch if (osm->lprolongation) PetscCall(VecScatterDestroy(&osm->lprolongation[i])); 6459566063dSJacob Faibussowitsch PetscCall(VecDestroy(&osm->x[i])); 6469566063dSJacob Faibussowitsch PetscCall(VecDestroy(&osm->y[i])); 6474b9ad928SBarry Smith } 6489566063dSJacob Faibussowitsch PetscCall(PetscFree(osm->lrestriction)); 6499566063dSJacob Faibussowitsch if (osm->lprolongation) PetscCall(PetscFree(osm->lprolongation)); 6509566063dSJacob Faibussowitsch PetscCall(PetscFree(osm->x)); 6519566063dSJacob Faibussowitsch PetscCall(PetscFree(osm->y)); 65292bb6962SLisandro Dalcin } 653ce78bad3SBarry Smith PetscCall(PCASMDestroySubdomains(osm->n_local_true, &osm->is, &osm->is_local)); 6549566063dSJacob Faibussowitsch PetscCall(ISDestroy(&osm->lis)); 6559566063dSJacob Faibussowitsch PetscCall(VecDestroy(&osm->lx)); 6569566063dSJacob Faibussowitsch PetscCall(VecDestroy(&osm->ly)); 65748a46eb9SPierre Jolivet if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE) PetscCall(MatDestroyMatrices(osm->n_local_true, &osm->lmats)); 6582fa5cd67SKarl Rupp 6599566063dSJacob Faibussowitsch PetscCall(PetscFree(osm->sub_mat_type)); 66080ec0b7dSPatrick Sanan 6610a545947SLisandro Dalcin osm->is = NULL; 6620a545947SLisandro Dalcin osm->is_local = NULL; 6633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 664e91c6855SBarry Smith } 665e91c6855SBarry Smith 666d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_ASM(PC pc) 667d71ae5a4SJacob Faibussowitsch { 668e91c6855SBarry Smith PC_ASM *osm = (PC_ASM *)pc->data; 669e91c6855SBarry Smith PetscInt i; 670e91c6855SBarry Smith 671e91c6855SBarry Smith PetscFunctionBegin; 6729566063dSJacob Faibussowitsch PetscCall(PCReset_ASM(pc)); 673e91c6855SBarry Smith if (osm->ksp) { 67448a46eb9SPierre Jolivet for (i = 0; i < osm->n_local_true; i++) PetscCall(KSPDestroy(&osm->ksp[i])); 6759566063dSJacob Faibussowitsch PetscCall(PetscFree(osm->ksp)); 676e91c6855SBarry Smith } 6779566063dSJacob Faibussowitsch PetscCall(PetscFree(pc->data)); 67896322394SPierre Jolivet 6799566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetLocalSubdomains_C", NULL)); 6809566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetTotalSubdomains_C", NULL)); 6819566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetOverlap_C", NULL)); 6829566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetType_C", NULL)); 6839566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMGetType_C", NULL)); 6849566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetLocalType_C", NULL)); 6859566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMGetLocalType_C", NULL)); 6869566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetSortIndices_C", NULL)); 6879566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMGetSubKSP_C", NULL)); 6889566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMGetSubMatType_C", NULL)); 6899566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetSubMatType_C", NULL)); 6903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6914b9ad928SBarry Smith } 6924b9ad928SBarry Smith 693ce78bad3SBarry Smith static PetscErrorCode PCSetFromOptions_ASM(PC pc, PetscOptionItems PetscOptionsObject) 694d71ae5a4SJacob Faibussowitsch { 6954b9ad928SBarry Smith PC_ASM *osm = (PC_ASM *)pc->data; 6969dcbbd2bSBarry Smith PetscInt blocks, ovl; 69757501b6eSBarry Smith PetscBool flg; 69892bb6962SLisandro Dalcin PCASMType asmtype; 69912cd4985SMatthew G. Knepley PCCompositeType loctype; 70080ec0b7dSPatrick Sanan char sub_mat_type[256]; 7014b9ad928SBarry Smith 7024b9ad928SBarry Smith PetscFunctionBegin; 703d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "Additive Schwarz options"); 7049566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_asm_dm_subdomains", "Use DMCreateDomainDecomposition() to define subdomains", "PCASMSetDMSubdomains", osm->dm_subdomains, &osm->dm_subdomains, &flg)); 7059566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_asm_blocks", "Number of subdomains", "PCASMSetTotalSubdomains", osm->n, &blocks, &flg)); 70665db9045SDmitry Karpeev if (flg) { 7079566063dSJacob Faibussowitsch PetscCall(PCASMSetTotalSubdomains(pc, blocks, NULL, NULL)); 708d709ea83SDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 70965db9045SDmitry Karpeev } 7109566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_asm_local_blocks", "Number of local subdomains", "PCASMSetLocalSubdomains", osm->n_local_true, &blocks, &flg)); 711342c94f9SMatthew G. Knepley if (flg) { 7129566063dSJacob Faibussowitsch PetscCall(PCASMSetLocalSubdomains(pc, blocks, NULL, NULL)); 713342c94f9SMatthew G. Knepley osm->dm_subdomains = PETSC_FALSE; 714342c94f9SMatthew G. Knepley } 7159566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_asm_overlap", "Number of grid points overlap", "PCASMSetOverlap", osm->overlap, &ovl, &flg)); 71665db9045SDmitry Karpeev if (flg) { 7179566063dSJacob Faibussowitsch PetscCall(PCASMSetOverlap(pc, ovl)); 718d709ea83SDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 71965db9045SDmitry Karpeev } 72090d69ab7SBarry Smith flg = PETSC_FALSE; 7219566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-pc_asm_type", "Type of restriction/extension", "PCASMSetType", PCASMTypes, (PetscEnum)osm->type, (PetscEnum *)&asmtype, &flg)); 7229566063dSJacob Faibussowitsch if (flg) PetscCall(PCASMSetType(pc, asmtype)); 72312cd4985SMatthew G. Knepley flg = PETSC_FALSE; 7249566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-pc_asm_local_type", "Type of local solver composition", "PCASMSetLocalType", PCCompositeTypes, (PetscEnum)osm->loctype, (PetscEnum *)&loctype, &flg)); 7259566063dSJacob Faibussowitsch if (flg) PetscCall(PCASMSetLocalType(pc, loctype)); 7269566063dSJacob Faibussowitsch PetscCall(PetscOptionsFList("-pc_asm_sub_mat_type", "Subsolve Matrix Type", "PCASMSetSubMatType", MatList, NULL, sub_mat_type, 256, &flg)); 7271baa6e33SBarry Smith if (flg) PetscCall(PCASMSetSubMatType(pc, sub_mat_type)); 728d0609cedSBarry Smith PetscOptionsHeadEnd(); 7293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7304b9ad928SBarry Smith } 7314b9ad928SBarry Smith 732d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCASMSetLocalSubdomains_ASM(PC pc, PetscInt n, IS is[], IS is_local[]) 733d71ae5a4SJacob Faibussowitsch { 7344b9ad928SBarry Smith PC_ASM *osm = (PC_ASM *)pc->data; 73592bb6962SLisandro Dalcin PetscInt i; 7364b9ad928SBarry Smith 7374b9ad928SBarry Smith PetscFunctionBegin; 73863a3b9bcSJacob Faibussowitsch PetscCheck(n >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Each process must have 1 or more blocks, n = %" PetscInt_FMT, n); 7397827d75bSBarry Smith PetscCheck(!pc->setupcalled || (n == osm->n_local_true && !is), PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "PCASMSetLocalSubdomains() should be called before calling PCSetUp()."); 740e7e72b3dSBarry Smith 7414b9ad928SBarry Smith if (!pc->setupcalled) { 74292bb6962SLisandro Dalcin if (is) { 7439566063dSJacob Faibussowitsch for (i = 0; i < n; i++) PetscCall(PetscObjectReference((PetscObject)is[i])); 74492bb6962SLisandro Dalcin } 745832fc9a5SMatthew Knepley if (is_local) { 7469566063dSJacob Faibussowitsch for (i = 0; i < n; i++) PetscCall(PetscObjectReference((PetscObject)is_local[i])); 747832fc9a5SMatthew Knepley } 748ce78bad3SBarry Smith PetscCall(PCASMDestroySubdomains(osm->n_local_true, &osm->is, &osm->is_local)); 7492fa5cd67SKarl Rupp 75007517c86SMark Adams if (osm->ksp && osm->n_local_true != n) { 75107517c86SMark Adams for (i = 0; i < osm->n_local_true; i++) PetscCall(KSPDestroy(&osm->ksp[i])); 75207517c86SMark Adams PetscCall(PetscFree(osm->ksp)); 75307517c86SMark Adams } 75407517c86SMark Adams 7554b9ad928SBarry Smith osm->n_local_true = n; 7560a545947SLisandro Dalcin osm->is = NULL; 7570a545947SLisandro Dalcin osm->is_local = NULL; 75892bb6962SLisandro Dalcin if (is) { 7599566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &osm->is)); 7602fa5cd67SKarl Rupp for (i = 0; i < n; i++) osm->is[i] = is[i]; 7613d03bb48SJed Brown /* Flag indicating that the user has set overlapping subdomains so PCASM should not increase their size. */ 7623d03bb48SJed Brown osm->overlap = -1; 76392bb6962SLisandro Dalcin } 7642b691e39SMatthew Knepley if (is_local) { 7659566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &osm->is_local)); 7662fa5cd67SKarl Rupp for (i = 0; i < n; i++) osm->is_local[i] = is_local[i]; 767a35b7b57SDmitry Karpeev if (!is) { 7689566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(osm->n_local_true, &osm->is)); 769a35b7b57SDmitry Karpeev for (i = 0; i < osm->n_local_true; i++) { 770a35b7b57SDmitry Karpeev if (osm->overlap > 0) { /* With positive overlap, osm->is[i] will be modified */ 7719566063dSJacob Faibussowitsch PetscCall(ISDuplicate(osm->is_local[i], &osm->is[i])); 7729566063dSJacob Faibussowitsch PetscCall(ISCopy(osm->is_local[i], osm->is[i])); 773a35b7b57SDmitry Karpeev } else { 7749566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)osm->is_local[i])); 775a35b7b57SDmitry Karpeev osm->is[i] = osm->is_local[i]; 776a35b7b57SDmitry Karpeev } 777a35b7b57SDmitry Karpeev } 778a35b7b57SDmitry Karpeev } 7792b691e39SMatthew Knepley } 7804b9ad928SBarry Smith } 7813ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7824b9ad928SBarry Smith } 7834b9ad928SBarry Smith 784d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCASMSetTotalSubdomains_ASM(PC pc, PetscInt N, IS *is, IS *is_local) 785d71ae5a4SJacob Faibussowitsch { 7864b9ad928SBarry Smith PC_ASM *osm = (PC_ASM *)pc->data; 78713f74950SBarry Smith PetscMPIInt rank, size; 78878904715SLisandro Dalcin PetscInt n; 7894b9ad928SBarry Smith 7904b9ad928SBarry Smith PetscFunctionBegin; 79163a3b9bcSJacob Faibussowitsch PetscCheck(N >= 1, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Number of total blocks must be > 0, N = %" PetscInt_FMT, N); 79200045ab3SPierre Jolivet PetscCheck(!is && !is_local, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Use PCASMSetLocalSubdomains() to set specific index sets, they cannot be set globally yet."); 7934b9ad928SBarry Smith 7944b9ad928SBarry Smith /* 795880770e9SJed Brown Split the subdomains equally among all processors 7964b9ad928SBarry Smith */ 7979566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank)); 7989566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size)); 7994b9ad928SBarry Smith n = N / size + ((N % size) > rank); 800835f2295SStefano 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); 8017827d75bSBarry Smith PetscCheck(!pc->setupcalled || n == osm->n_local_true, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "PCASMSetTotalSubdomains() should be called before PCSetUp()."); 8024b9ad928SBarry Smith if (!pc->setupcalled) { 803ce78bad3SBarry Smith PetscCall(PCASMDestroySubdomains(osm->n_local_true, &osm->is, &osm->is_local)); 8042fa5cd67SKarl Rupp 8054b9ad928SBarry Smith osm->n_local_true = n; 8060a545947SLisandro Dalcin osm->is = NULL; 8070a545947SLisandro Dalcin osm->is_local = NULL; 8084b9ad928SBarry Smith } 8093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8104b9ad928SBarry Smith } 8114b9ad928SBarry Smith 812d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCASMSetOverlap_ASM(PC pc, PetscInt ovl) 813d71ae5a4SJacob Faibussowitsch { 81492bb6962SLisandro Dalcin PC_ASM *osm = (PC_ASM *)pc->data; 8154b9ad928SBarry Smith 8164b9ad928SBarry Smith PetscFunctionBegin; 8177827d75bSBarry Smith PetscCheck(ovl >= 0, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Negative overlap value requested"); 8187827d75bSBarry Smith PetscCheck(!pc->setupcalled || ovl == osm->overlap, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "PCASMSetOverlap() should be called before PCSetUp()."); 8192fa5cd67SKarl Rupp if (!pc->setupcalled) osm->overlap = ovl; 8203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8214b9ad928SBarry Smith } 8224b9ad928SBarry Smith 823d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCASMSetType_ASM(PC pc, PCASMType type) 824d71ae5a4SJacob Faibussowitsch { 82592bb6962SLisandro Dalcin PC_ASM *osm = (PC_ASM *)pc->data; 8264b9ad928SBarry Smith 8274b9ad928SBarry Smith PetscFunctionBegin; 8284b9ad928SBarry Smith osm->type = type; 829bf108f30SBarry Smith osm->type_set = PETSC_TRUE; 8303ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8314b9ad928SBarry Smith } 8324b9ad928SBarry Smith 833d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCASMGetType_ASM(PC pc, PCASMType *type) 834d71ae5a4SJacob Faibussowitsch { 835c60c7ad4SBarry Smith PC_ASM *osm = (PC_ASM *)pc->data; 836c60c7ad4SBarry Smith 837c60c7ad4SBarry Smith PetscFunctionBegin; 838c60c7ad4SBarry Smith *type = osm->type; 8393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 840c60c7ad4SBarry Smith } 841c60c7ad4SBarry Smith 842d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCASMSetLocalType_ASM(PC pc, PCCompositeType type) 843d71ae5a4SJacob Faibussowitsch { 84412cd4985SMatthew G. Knepley PC_ASM *osm = (PC_ASM *)pc->data; 84512cd4985SMatthew G. Knepley 84612cd4985SMatthew G. Knepley PetscFunctionBegin; 8477827d75bSBarry 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"); 84812cd4985SMatthew G. Knepley osm->loctype = type; 8493ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 85012cd4985SMatthew G. Knepley } 85112cd4985SMatthew G. Knepley 852d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCASMGetLocalType_ASM(PC pc, PCCompositeType *type) 853d71ae5a4SJacob Faibussowitsch { 85412cd4985SMatthew G. Knepley PC_ASM *osm = (PC_ASM *)pc->data; 85512cd4985SMatthew G. Knepley 85612cd4985SMatthew G. Knepley PetscFunctionBegin; 85712cd4985SMatthew G. Knepley *type = osm->loctype; 8583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 85912cd4985SMatthew G. Knepley } 86012cd4985SMatthew G. Knepley 861d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCASMSetSortIndices_ASM(PC pc, PetscBool doSort) 862d71ae5a4SJacob Faibussowitsch { 8636ed231c7SMatthew Knepley PC_ASM *osm = (PC_ASM *)pc->data; 8646ed231c7SMatthew Knepley 8656ed231c7SMatthew Knepley PetscFunctionBegin; 8666ed231c7SMatthew Knepley osm->sort_indices = doSort; 8673ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8686ed231c7SMatthew Knepley } 8696ed231c7SMatthew Knepley 870d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCASMGetSubKSP_ASM(PC pc, PetscInt *n_local, PetscInt *first_local, KSP **ksp) 871d71ae5a4SJacob Faibussowitsch { 87292bb6962SLisandro Dalcin PC_ASM *osm = (PC_ASM *)pc->data; 8734b9ad928SBarry Smith 8744b9ad928SBarry Smith PetscFunctionBegin; 8757827d75bSBarry 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"); 8764b9ad928SBarry Smith 8772fa5cd67SKarl Rupp if (n_local) *n_local = osm->n_local_true; 87892bb6962SLisandro Dalcin if (first_local) { 8799566063dSJacob Faibussowitsch PetscCallMPI(MPI_Scan(&osm->n_local_true, first_local, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)pc))); 88092bb6962SLisandro Dalcin *first_local -= osm->n_local_true; 88192bb6962SLisandro Dalcin } 882ed155784SPierre Jolivet if (ksp) *ksp = osm->ksp; 8833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8844b9ad928SBarry Smith } 8854b9ad928SBarry Smith 886d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCASMGetSubMatType_ASM(PC pc, MatType *sub_mat_type) 887d71ae5a4SJacob Faibussowitsch { 88880ec0b7dSPatrick Sanan PC_ASM *osm = (PC_ASM *)pc->data; 88980ec0b7dSPatrick Sanan 89080ec0b7dSPatrick Sanan PetscFunctionBegin; 89180ec0b7dSPatrick Sanan PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 8924f572ea9SToby Isaac PetscAssertPointer(sub_mat_type, 2); 89380ec0b7dSPatrick Sanan *sub_mat_type = osm->sub_mat_type; 8943ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 89580ec0b7dSPatrick Sanan } 89680ec0b7dSPatrick Sanan 897d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCASMSetSubMatType_ASM(PC pc, MatType sub_mat_type) 898d71ae5a4SJacob Faibussowitsch { 89980ec0b7dSPatrick Sanan PC_ASM *osm = (PC_ASM *)pc->data; 90080ec0b7dSPatrick Sanan 90180ec0b7dSPatrick Sanan PetscFunctionBegin; 90280ec0b7dSPatrick Sanan PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 9039566063dSJacob Faibussowitsch PetscCall(PetscFree(osm->sub_mat_type)); 9049566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(sub_mat_type, (char **)&osm->sub_mat_type)); 9053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 90680ec0b7dSPatrick Sanan } 90780ec0b7dSPatrick Sanan 9085d83a8b1SBarry Smith /*@ 909f1580f4eSBarry Smith PCASMSetLocalSubdomains - Sets the local subdomains (for this processor only) for the additive Schwarz preconditioner `PCASM`. 9104b9ad928SBarry Smith 911c3339decSBarry Smith Collective 9124b9ad928SBarry Smith 9134b9ad928SBarry Smith Input Parameters: 9144b9ad928SBarry Smith + pc - the preconditioner context 9154b9ad928SBarry Smith . n - the number of subdomains for this processor (default value = 1) 916a3b724e8SBarry Smith . is - the index set that defines the subdomains for this processor (or `NULL` for PETSc to determine subdomains) 917f13dfd9eSBarry Smith the values of the `is` array are copied so you can free the array (not the `IS` in the array) after this call 918f13dfd9eSBarry 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` 919f13dfd9eSBarry Smith (or `NULL` to not provide these). The values of the `is_local` array are copied so you can free the array 920f13dfd9eSBarry Smith (not the `IS` in the array) after this call 9214b9ad928SBarry Smith 922342c94f9SMatthew G. Knepley Options Database Key: 92320f4b53cSBarry Smith . -pc_asm_local_blocks <blks> - Sets number of local blocks 92420f4b53cSBarry Smith 92520f4b53cSBarry Smith Level: advanced 926342c94f9SMatthew G. Knepley 9274b9ad928SBarry Smith Notes: 928a3b724e8SBarry Smith The `IS` numbering is in the parallel, global numbering of the vector for both `is` and `is_local` 9294b9ad928SBarry Smith 930f1580f4eSBarry Smith By default the `PCASM` preconditioner uses 1 block per processor. 9314b9ad928SBarry Smith 932f1580f4eSBarry Smith Use `PCASMSetTotalSubdomains()` to set the subdomains for all processors. 9334b9ad928SBarry Smith 934a3b724e8SBarry Smith If `is_local` is provided and `PCASMType` is `PC_ASM_RESTRICT` then the solution only over the `is_local` region is interpolated 935a3b724e8SBarry Smith back to form the global solution (this is the standard restricted additive Schwarz method, RASM) 936f1ee410cSBarry Smith 937a3b724e8SBarry Smith If `is_local` is provided and `PCASMType` is `PC_ASM_INTERPOLATE` or `PC_ASM_NONE` then an error is generated since there is 938f1ee410cSBarry Smith no code to handle that case. 939f1ee410cSBarry Smith 940562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMSetTotalSubdomains()`, `PCASMSetOverlap()`, `PCASMGetSubKSP()`, 941f1580f4eSBarry Smith `PCASMCreateSubdomains2D()`, `PCASMGetLocalSubdomains()`, `PCASMType`, `PCASMSetType()`, `PCGASM` 9424b9ad928SBarry Smith @*/ 943d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMSetLocalSubdomains(PC pc, PetscInt n, IS is[], IS is_local[]) 944d71ae5a4SJacob Faibussowitsch { 9454b9ad928SBarry Smith PetscFunctionBegin; 9460700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 947cac4c232SBarry Smith PetscTryMethod(pc, "PCASMSetLocalSubdomains_C", (PC, PetscInt, IS[], IS[]), (pc, n, is, is_local)); 9483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9494b9ad928SBarry Smith } 9504b9ad928SBarry Smith 9515d83a8b1SBarry Smith /*@ 952feb221f8SDmitry Karpeev PCASMSetTotalSubdomains - Sets the subdomains for all processors for the 953f1580f4eSBarry Smith additive Schwarz preconditioner, `PCASM`. 9544b9ad928SBarry Smith 955c3339decSBarry Smith Collective, all MPI ranks must pass in the same array of `IS` 9564b9ad928SBarry Smith 9574b9ad928SBarry Smith Input Parameters: 9584b9ad928SBarry Smith + pc - the preconditioner context 959feb221f8SDmitry Karpeev . N - the number of subdomains for all processors 960a3b724e8SBarry Smith . is - the index sets that define the subdomains for all processors (or `NULL` to ask PETSc to determine the subdomains) 9615d83a8b1SBarry Smith the values of the `is` array are copied so you can free the array (not the `IS` in the array) after this call 962a3b724e8SBarry Smith - is_local - the index sets that define the local part of the subdomains for this processor (or `NULL` to not provide this information) 9635d83a8b1SBarry 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 9644b9ad928SBarry Smith 9654b9ad928SBarry Smith Options Database Key: 9664b9ad928SBarry Smith . -pc_asm_blocks <blks> - Sets total blocks 9674b9ad928SBarry Smith 96820f4b53cSBarry Smith Level: advanced 96920f4b53cSBarry Smith 9704b9ad928SBarry Smith Notes: 971a3b724e8SBarry Smith Currently you cannot use this to set the actual subdomains with the argument `is` or `is_local`. 9724b9ad928SBarry Smith 973f1580f4eSBarry Smith By default the `PCASM` preconditioner uses 1 block per processor. 9744b9ad928SBarry Smith 9754b9ad928SBarry Smith These index sets cannot be destroyed until after completion of the 976f1580f4eSBarry Smith linear solves for which the `PCASM` preconditioner is being used. 9774b9ad928SBarry Smith 978f1580f4eSBarry Smith Use `PCASMSetLocalSubdomains()` to set local subdomains. 9794b9ad928SBarry Smith 980f1580f4eSBarry Smith The `IS` numbering is in the parallel, global numbering of the vector for both is and is_local 9811093a601SBarry Smith 982562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMSetLocalSubdomains()`, `PCASMSetOverlap()`, `PCASMGetSubKSP()`, 983f1580f4eSBarry Smith `PCASMCreateSubdomains2D()`, `PCGASM` 9844b9ad928SBarry Smith @*/ 985d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMSetTotalSubdomains(PC pc, PetscInt N, IS is[], IS is_local[]) 986d71ae5a4SJacob Faibussowitsch { 9874b9ad928SBarry Smith PetscFunctionBegin; 9880700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 989cac4c232SBarry Smith PetscTryMethod(pc, "PCASMSetTotalSubdomains_C", (PC, PetscInt, IS[], IS[]), (pc, N, is, is_local)); 9903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9914b9ad928SBarry Smith } 9924b9ad928SBarry Smith 9934b9ad928SBarry Smith /*@ 9944b9ad928SBarry Smith PCASMSetOverlap - Sets the overlap between a pair of subdomains for the 995f1580f4eSBarry Smith additive Schwarz preconditioner, `PCASM`. 9964b9ad928SBarry Smith 997c3339decSBarry Smith Logically Collective 9984b9ad928SBarry Smith 9994b9ad928SBarry Smith Input Parameters: 10004b9ad928SBarry Smith + pc - the preconditioner context 10014b9ad928SBarry Smith - ovl - the amount of overlap between subdomains (ovl >= 0, default value = 1) 10024b9ad928SBarry Smith 10034b9ad928SBarry Smith Options Database Key: 10044b9ad928SBarry Smith . -pc_asm_overlap <ovl> - Sets overlap 10054b9ad928SBarry Smith 100620f4b53cSBarry Smith Level: intermediate 100720f4b53cSBarry Smith 10084b9ad928SBarry Smith Notes: 1009f1580f4eSBarry Smith By default the `PCASM` preconditioner uses 1 block per processor. To use 1010f1580f4eSBarry Smith multiple blocks per perocessor, see `PCASMSetTotalSubdomains()` and 1011f1580f4eSBarry Smith `PCASMSetLocalSubdomains()` (and the option -pc_asm_blocks <blks>). 10124b9ad928SBarry Smith 10134b9ad928SBarry Smith The overlap defaults to 1, so if one desires that no additional 10144b9ad928SBarry Smith overlap be computed beyond what may have been set with a call to 1015f1580f4eSBarry Smith `PCASMSetTotalSubdomains()` or `PCASMSetLocalSubdomains()`, then ovl 10164b9ad928SBarry Smith must be set to be 0. In particular, if one does not explicitly set 10174b9ad928SBarry Smith the subdomains an application code, then all overlap would be computed 1018f1580f4eSBarry Smith internally by PETSc, and using an overlap of 0 would result in an `PCASM` 10194b9ad928SBarry Smith variant that is equivalent to the block Jacobi preconditioner. 10204b9ad928SBarry Smith 1021f1ee410cSBarry Smith The default algorithm used by PETSc to increase overlap is fast, but not scalable, 1022f1ee410cSBarry Smith use the option -mat_increase_overlap_scalable when the problem and number of processes is large. 1023f1ee410cSBarry Smith 10242fe279fdSBarry Smith One can define initial index sets with any overlap via 1025f1580f4eSBarry Smith `PCASMSetLocalSubdomains()`; the routine 1026f1580f4eSBarry Smith `PCASMSetOverlap()` merely allows PETSc to extend that overlap further 10274b9ad928SBarry Smith if desired. 10284b9ad928SBarry Smith 1029562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMSetTotalSubdomains()`, `PCASMSetLocalSubdomains()`, `PCASMGetSubKSP()`, 1030f1580f4eSBarry Smith `PCASMCreateSubdomains2D()`, `PCASMGetLocalSubdomains()`, `MatIncreaseOverlap()`, `PCGASM` 10314b9ad928SBarry Smith @*/ 1032d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMSetOverlap(PC pc, PetscInt ovl) 1033d71ae5a4SJacob Faibussowitsch { 10344b9ad928SBarry Smith PetscFunctionBegin; 10350700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1036c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc, ovl, 2); 1037cac4c232SBarry Smith PetscTryMethod(pc, "PCASMSetOverlap_C", (PC, PetscInt), (pc, ovl)); 10383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 10394b9ad928SBarry Smith } 10404b9ad928SBarry Smith 10414b9ad928SBarry Smith /*@ 10424b9ad928SBarry Smith PCASMSetType - Sets the type of restriction and interpolation used 1043f1580f4eSBarry Smith for local problems in the additive Schwarz method, `PCASM`. 10444b9ad928SBarry Smith 1045c3339decSBarry Smith Logically Collective 10464b9ad928SBarry Smith 10474b9ad928SBarry Smith Input Parameters: 10484b9ad928SBarry Smith + pc - the preconditioner context 1049f1580f4eSBarry Smith - type - variant of `PCASM`, one of 10504b9ad928SBarry Smith .vb 10514b9ad928SBarry Smith PC_ASM_BASIC - full interpolation and restriction 105282b5ce2aSStefano Zampini PC_ASM_RESTRICT - full restriction, local processor interpolation (default) 10534b9ad928SBarry Smith PC_ASM_INTERPOLATE - full interpolation, local processor restriction 10544b9ad928SBarry Smith PC_ASM_NONE - local processor restriction and interpolation 10554b9ad928SBarry Smith .ve 10564b9ad928SBarry Smith 10574b9ad928SBarry Smith Options Database Key: 1058f1580f4eSBarry Smith . -pc_asm_type [basic,restrict,interpolate,none] - Sets `PCASMType` 10594b9ad928SBarry Smith 106020f4b53cSBarry Smith Level: intermediate 106120f4b53cSBarry Smith 1062f1580f4eSBarry Smith Note: 1063f1580f4eSBarry Smith if the is_local arguments are passed to `PCASMSetLocalSubdomains()` then they are used when `PC_ASM_RESTRICT` has been selected 1064f1ee410cSBarry Smith to limit the local processor interpolation 1065f1ee410cSBarry Smith 1066562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMSetTotalSubdomains()`, `PCASMGetSubKSP()`, 1067f1580f4eSBarry Smith `PCASMCreateSubdomains2D()`, `PCASMType`, `PCASMSetLocalType()`, `PCASMGetLocalType()`, `PCGASM` 10684b9ad928SBarry Smith @*/ 1069d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMSetType(PC pc, PCASMType type) 1070d71ae5a4SJacob Faibussowitsch { 10714b9ad928SBarry Smith PetscFunctionBegin; 10720700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1073c5eb9154SBarry Smith PetscValidLogicalCollectiveEnum(pc, type, 2); 1074cac4c232SBarry Smith PetscTryMethod(pc, "PCASMSetType_C", (PC, PCASMType), (pc, type)); 10753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 10764b9ad928SBarry Smith } 10774b9ad928SBarry Smith 1078c60c7ad4SBarry Smith /*@ 1079c60c7ad4SBarry Smith PCASMGetType - Gets the type of restriction and interpolation used 1080f1580f4eSBarry Smith for local problems in the additive Schwarz method, `PCASM`. 1081c60c7ad4SBarry Smith 1082c3339decSBarry Smith Logically Collective 1083c60c7ad4SBarry Smith 1084c60c7ad4SBarry Smith Input Parameter: 1085c60c7ad4SBarry Smith . pc - the preconditioner context 1086c60c7ad4SBarry Smith 1087c60c7ad4SBarry Smith Output Parameter: 1088f1580f4eSBarry Smith . type - variant of `PCASM`, one of 1089c60c7ad4SBarry Smith .vb 1090c60c7ad4SBarry Smith PC_ASM_BASIC - full interpolation and restriction 1091c60c7ad4SBarry Smith PC_ASM_RESTRICT - full restriction, local processor interpolation 1092c60c7ad4SBarry Smith PC_ASM_INTERPOLATE - full interpolation, local processor restriction 1093c60c7ad4SBarry Smith PC_ASM_NONE - local processor restriction and interpolation 1094c60c7ad4SBarry Smith .ve 1095c60c7ad4SBarry Smith 1096c60c7ad4SBarry Smith Options Database Key: 1097f1580f4eSBarry Smith . -pc_asm_type [basic,restrict,interpolate,none] - Sets `PCASM` type 1098c60c7ad4SBarry Smith 1099c60c7ad4SBarry Smith Level: intermediate 1100c60c7ad4SBarry Smith 1101562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMSetTotalSubdomains()`, `PCASMGetSubKSP()`, `PCGASM`, 1102db781477SPatrick Sanan `PCASMCreateSubdomains2D()`, `PCASMType`, `PCASMSetType()`, `PCASMSetLocalType()`, `PCASMGetLocalType()` 1103c60c7ad4SBarry Smith @*/ 1104d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMGetType(PC pc, PCASMType *type) 1105d71ae5a4SJacob Faibussowitsch { 1106c60c7ad4SBarry Smith PetscFunctionBegin; 1107c60c7ad4SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1108cac4c232SBarry Smith PetscUseMethod(pc, "PCASMGetType_C", (PC, PCASMType *), (pc, type)); 11093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1110c60c7ad4SBarry Smith } 1111c60c7ad4SBarry Smith 111212cd4985SMatthew G. Knepley /*@ 1113f1580f4eSBarry Smith PCASMSetLocalType - Sets the type of composition used for local problems in the additive Schwarz method, `PCASM`. 111412cd4985SMatthew G. Knepley 1115c3339decSBarry Smith Logically Collective 111612cd4985SMatthew G. Knepley 111712cd4985SMatthew G. Knepley Input Parameters: 111812cd4985SMatthew G. Knepley + pc - the preconditioner context 111912cd4985SMatthew G. Knepley - type - type of composition, one of 112012cd4985SMatthew G. Knepley .vb 112112cd4985SMatthew G. Knepley PC_COMPOSITE_ADDITIVE - local additive combination 112212cd4985SMatthew G. Knepley PC_COMPOSITE_MULTIPLICATIVE - local multiplicative combination 112312cd4985SMatthew G. Knepley .ve 112412cd4985SMatthew G. Knepley 112512cd4985SMatthew G. Knepley Options Database Key: 112612cd4985SMatthew G. Knepley . -pc_asm_local_type [additive,multiplicative] - Sets local solver composition type 112712cd4985SMatthew G. Knepley 112812cd4985SMatthew G. Knepley Level: intermediate 112912cd4985SMatthew G. Knepley 1130562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMSetType()`, `PCASMGetType()`, `PCASMGetLocalType()`, `PCASMType`, `PCCompositeType` 113112cd4985SMatthew G. Knepley @*/ 1132d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMSetLocalType(PC pc, PCCompositeType type) 1133d71ae5a4SJacob Faibussowitsch { 113412cd4985SMatthew G. Knepley PetscFunctionBegin; 113512cd4985SMatthew G. Knepley PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 113612cd4985SMatthew G. Knepley PetscValidLogicalCollectiveEnum(pc, type, 2); 1137cac4c232SBarry Smith PetscTryMethod(pc, "PCASMSetLocalType_C", (PC, PCCompositeType), (pc, type)); 11383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 113912cd4985SMatthew G. Knepley } 114012cd4985SMatthew G. Knepley 114112cd4985SMatthew G. Knepley /*@ 1142f1580f4eSBarry Smith PCASMGetLocalType - Gets the type of composition used for local problems in the additive Schwarz method, `PCASM`. 114312cd4985SMatthew G. Knepley 1144c3339decSBarry Smith Logically Collective 114512cd4985SMatthew G. Knepley 114612cd4985SMatthew G. Knepley Input Parameter: 114712cd4985SMatthew G. Knepley . pc - the preconditioner context 114812cd4985SMatthew G. Knepley 114912cd4985SMatthew G. Knepley Output Parameter: 115012cd4985SMatthew G. Knepley . type - type of composition, one of 115112cd4985SMatthew G. Knepley .vb 115212cd4985SMatthew G. Knepley PC_COMPOSITE_ADDITIVE - local additive combination 115312cd4985SMatthew G. Knepley PC_COMPOSITE_MULTIPLICATIVE - local multiplicative combination 115412cd4985SMatthew G. Knepley .ve 115512cd4985SMatthew G. Knepley 115612cd4985SMatthew G. Knepley Options Database Key: 115712cd4985SMatthew G. Knepley . -pc_asm_local_type [additive,multiplicative] - Sets local solver composition type 115812cd4985SMatthew G. Knepley 115912cd4985SMatthew G. Knepley Level: intermediate 116012cd4985SMatthew G. Knepley 11610241b274SPierre Jolivet .seealso: [](ch_ksp), `PCASM`, `PCASMSetType()`, `PCASMGetType()`, `PCASMSetLocalType()`, `PCASMType`, `PCCompositeType` 116212cd4985SMatthew G. Knepley @*/ 1163d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMGetLocalType(PC pc, PCCompositeType *type) 1164d71ae5a4SJacob Faibussowitsch { 116512cd4985SMatthew G. Knepley PetscFunctionBegin; 116612cd4985SMatthew G. Knepley PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 11674f572ea9SToby Isaac PetscAssertPointer(type, 2); 1168cac4c232SBarry Smith PetscUseMethod(pc, "PCASMGetLocalType_C", (PC, PCCompositeType *), (pc, type)); 11693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 117012cd4985SMatthew G. Knepley } 117112cd4985SMatthew G. Knepley 11726ed231c7SMatthew Knepley /*@ 11736ed231c7SMatthew Knepley PCASMSetSortIndices - Determines whether subdomain indices are sorted. 11746ed231c7SMatthew Knepley 1175c3339decSBarry Smith Logically Collective 11766ed231c7SMatthew Knepley 11776ed231c7SMatthew Knepley Input Parameters: 11786ed231c7SMatthew Knepley + pc - the preconditioner context 11796ed231c7SMatthew Knepley - doSort - sort the subdomain indices 11806ed231c7SMatthew Knepley 11816ed231c7SMatthew Knepley Level: intermediate 11826ed231c7SMatthew Knepley 1183562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMSetLocalSubdomains()`, `PCASMSetTotalSubdomains()`, `PCASMGetSubKSP()`, 1184db781477SPatrick Sanan `PCASMCreateSubdomains2D()` 11856ed231c7SMatthew Knepley @*/ 1186d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMSetSortIndices(PC pc, PetscBool doSort) 1187d71ae5a4SJacob Faibussowitsch { 11886ed231c7SMatthew Knepley PetscFunctionBegin; 11890700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1190acfcf0e5SJed Brown PetscValidLogicalCollectiveBool(pc, doSort, 2); 1191cac4c232SBarry Smith PetscTryMethod(pc, "PCASMSetSortIndices_C", (PC, PetscBool), (pc, doSort)); 11923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11936ed231c7SMatthew Knepley } 11946ed231c7SMatthew Knepley 11954b9ad928SBarry Smith /*@C 1196f1580f4eSBarry Smith PCASMGetSubKSP - Gets the local `KSP` contexts for all blocks on 11974b9ad928SBarry Smith this processor. 11984b9ad928SBarry Smith 1199c3339decSBarry Smith Collective iff first_local is requested 12004b9ad928SBarry Smith 12014b9ad928SBarry Smith Input Parameter: 12024b9ad928SBarry Smith . pc - the preconditioner context 12034b9ad928SBarry Smith 12044b9ad928SBarry Smith Output Parameters: 1205a3b724e8SBarry Smith + n_local - the number of blocks on this processor or `NULL` 1206a3b724e8SBarry Smith . first_local - the global number of the first block on this processor or `NULL`, all processors must request or all must pass `NULL` 1207f1580f4eSBarry Smith - ksp - the array of `KSP` contexts 12084b9ad928SBarry Smith 120920f4b53cSBarry Smith Level: advanced 121020f4b53cSBarry Smith 1211f1580f4eSBarry Smith Notes: 1212f1580f4eSBarry Smith After `PCASMGetSubKSP()` the array of `KSP`s is not to be freed. 12134b9ad928SBarry Smith 1214f1580f4eSBarry Smith You must call `KSPSetUp()` before calling `PCASMGetSubKSP()`. 12154b9ad928SBarry Smith 1216562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMSetTotalSubdomains()`, `PCASMSetOverlap()`, 1217db781477SPatrick Sanan `PCASMCreateSubdomains2D()`, 12184b9ad928SBarry Smith @*/ 1219d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMGetSubKSP(PC pc, PetscInt *n_local, PetscInt *first_local, KSP *ksp[]) 1220d71ae5a4SJacob Faibussowitsch { 12214b9ad928SBarry Smith PetscFunctionBegin; 12220700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1223cac4c232SBarry Smith PetscUseMethod(pc, "PCASMGetSubKSP_C", (PC, PetscInt *, PetscInt *, KSP **), (pc, n_local, first_local, ksp)); 12243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 12254b9ad928SBarry Smith } 12264b9ad928SBarry Smith 12274b9ad928SBarry Smith /*MC 12283b09bd56SBarry Smith PCASM - Use the (restricted) additive Schwarz method, each block is (approximately) solved with 12291d27aa22SBarry Smith its own `KSP` object, {cite}`dryja1987additive` and {cite}`1sbg` 12304b9ad928SBarry Smith 12314b9ad928SBarry Smith Options Database Keys: 12321d27aa22SBarry Smith + -pc_asm_blocks <blks> - Sets total blocks. Defaults to one block per MPI process. 12334b9ad928SBarry Smith . -pc_asm_overlap <ovl> - Sets overlap 1234f1580f4eSBarry Smith . -pc_asm_type [basic,restrict,interpolate,none] - Sets `PCASMType`, default is restrict. See `PCASMSetType()` 123573ff1848SBarry Smith . -pc_asm_dm_subdomains <bool> - use subdomains defined by the `DM` with `DMCreateDomainDecomposition()` 1236f1580f4eSBarry Smith - -pc_asm_local_type [additive, multiplicative] - Sets `PCCompositeType`, default is additive. See `PCASMSetLocalType()` 12374b9ad928SBarry Smith 12384b9ad928SBarry Smith Level: beginner 12394b9ad928SBarry Smith 1240f1580f4eSBarry Smith Notes: 1241f1580f4eSBarry Smith If you run with, for example, 3 blocks on 1 processor or 3 blocks on 3 processors you 124273ff1848SBarry Smith will get a different convergence rate due to the default option of `-pc_asm_type restrict`. Use 124373ff1848SBarry Smith `-pc_asm_type basic` to get the same convergence behavior 1244f1580f4eSBarry Smith 1245f1580f4eSBarry Smith Each processor can have one or more blocks, but a block cannot be shared by more 1246f1580f4eSBarry Smith than one processor. Use `PCGASM` for subdomains shared by multiple processes. 1247f1580f4eSBarry Smith 124873ff1848SBarry Smith To set options on the solvers for each block append `-sub_` to all the `KSP`, and `PC` 124973ff1848SBarry Smith options database keys. For example, `-sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly` 1250f1580f4eSBarry Smith 1251f1580f4eSBarry Smith To set the options on the solvers separate for each block call `PCASMGetSubKSP()` 1252f1580f4eSBarry Smith and set the options directly on the resulting `KSP` object (you can access its `PC` with `KSPGetPC()`) 1253f1580f4eSBarry Smith 125473ff1848SBarry Smith If the `PC` has an associated `DM`, then, by default, `DMCreateDomainDecomposition()` is used to create the subdomains 125573ff1848SBarry Smith 1256562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCASMType`, `PCCompositeType`, 1257db781477SPatrick Sanan `PCBJACOBI`, `PCASMGetSubKSP()`, `PCASMSetLocalSubdomains()`, `PCASMType`, `PCASMGetType()`, `PCASMSetLocalType()`, `PCASMGetLocalType()` 1258db781477SPatrick Sanan `PCASMSetTotalSubdomains()`, `PCSetModifySubMatrices()`, `PCASMSetOverlap()`, `PCASMSetType()`, `PCCompositeType` 12594b9ad928SBarry Smith M*/ 12604b9ad928SBarry Smith 1261d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_ASM(PC pc) 1262d71ae5a4SJacob Faibussowitsch { 12634b9ad928SBarry Smith PC_ASM *osm; 12644b9ad928SBarry Smith 12654b9ad928SBarry Smith PetscFunctionBegin; 12664dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&osm)); 12672fa5cd67SKarl Rupp 12684b9ad928SBarry Smith osm->n = PETSC_DECIDE; 12694b9ad928SBarry Smith osm->n_local = 0; 12702b837212SDmitry Karpeev osm->n_local_true = PETSC_DECIDE; 12714b9ad928SBarry Smith osm->overlap = 1; 12720a545947SLisandro Dalcin osm->ksp = NULL; 12730a545947SLisandro Dalcin osm->restriction = NULL; 12740a545947SLisandro Dalcin osm->lprolongation = NULL; 12750a545947SLisandro Dalcin osm->lrestriction = NULL; 12760a545947SLisandro Dalcin osm->x = NULL; 12770a545947SLisandro Dalcin osm->y = NULL; 12780a545947SLisandro Dalcin osm->is = NULL; 12790a545947SLisandro Dalcin osm->is_local = NULL; 12800a545947SLisandro Dalcin osm->mat = NULL; 12810a545947SLisandro Dalcin osm->pmat = NULL; 12824b9ad928SBarry Smith osm->type = PC_ASM_RESTRICT; 128312cd4985SMatthew G. Knepley osm->loctype = PC_COMPOSITE_ADDITIVE; 12846ed231c7SMatthew Knepley osm->sort_indices = PETSC_TRUE; 1285d709ea83SDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 128680ec0b7dSPatrick Sanan osm->sub_mat_type = NULL; 12874b9ad928SBarry Smith 128892bb6962SLisandro Dalcin pc->data = (void *)osm; 12894b9ad928SBarry Smith pc->ops->apply = PCApply_ASM; 129048e38a8aSPierre Jolivet pc->ops->matapply = PCMatApply_ASM; 12914b9ad928SBarry Smith pc->ops->applytranspose = PCApplyTranspose_ASM; 12924b9ad928SBarry Smith pc->ops->setup = PCSetUp_ASM; 1293e91c6855SBarry Smith pc->ops->reset = PCReset_ASM; 12944b9ad928SBarry Smith pc->ops->destroy = PCDestroy_ASM; 12954b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_ASM; 12964b9ad928SBarry Smith pc->ops->setuponblocks = PCSetUpOnBlocks_ASM; 12974b9ad928SBarry Smith pc->ops->view = PCView_ASM; 12980a545947SLisandro Dalcin pc->ops->applyrichardson = NULL; 12994b9ad928SBarry Smith 13009566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetLocalSubdomains_C", PCASMSetLocalSubdomains_ASM)); 13019566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetTotalSubdomains_C", PCASMSetTotalSubdomains_ASM)); 13029566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetOverlap_C", PCASMSetOverlap_ASM)); 13039566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetType_C", PCASMSetType_ASM)); 13049566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMGetType_C", PCASMGetType_ASM)); 13059566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetLocalType_C", PCASMSetLocalType_ASM)); 13069566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMGetLocalType_C", PCASMGetLocalType_ASM)); 13079566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetSortIndices_C", PCASMSetSortIndices_ASM)); 13089566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMGetSubKSP_C", PCASMGetSubKSP_ASM)); 13099566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMGetSubMatType_C", PCASMGetSubMatType_ASM)); 13109566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetSubMatType_C", PCASMSetSubMatType_ASM)); 13113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13124b9ad928SBarry Smith } 13134b9ad928SBarry Smith 131492bb6962SLisandro Dalcin /*@C 131592bb6962SLisandro Dalcin PCASMCreateSubdomains - Creates the index sets for the overlapping Schwarz 1316f1580f4eSBarry Smith preconditioner, `PCASM`, for any problem on a general grid. 131792bb6962SLisandro Dalcin 131892bb6962SLisandro Dalcin Collective 131992bb6962SLisandro Dalcin 132092bb6962SLisandro Dalcin Input Parameters: 132192bb6962SLisandro Dalcin + A - The global matrix operator 132292bb6962SLisandro Dalcin - n - the number of local blocks 132392bb6962SLisandro Dalcin 13242fe279fdSBarry Smith Output Parameter: 132592bb6962SLisandro Dalcin . outis - the array of index sets defining the subdomains 132692bb6962SLisandro Dalcin 132792bb6962SLisandro Dalcin Level: advanced 132892bb6962SLisandro Dalcin 1329f1580f4eSBarry Smith Note: 1330f1580f4eSBarry Smith This generates nonoverlapping subdomains; the `PCASM` will generate the overlap 1331f1580f4eSBarry Smith from these if you use `PCASMSetLocalSubdomains()` 13327d6bfa3bSBarry Smith 1333562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMSetLocalSubdomains()`, `PCASMDestroySubdomains()` 133492bb6962SLisandro Dalcin @*/ 1335d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMCreateSubdomains(Mat A, PetscInt n, IS *outis[]) 1336d71ae5a4SJacob Faibussowitsch { 133792bb6962SLisandro Dalcin MatPartitioning mpart; 133892bb6962SLisandro Dalcin const char *prefix; 133992bb6962SLisandro Dalcin PetscInt i, j, rstart, rend, bs; 1340976e8c5aSLisandro Dalcin PetscBool hasop, isbaij = PETSC_FALSE, foundpart = PETSC_FALSE; 13410298fd71SBarry Smith Mat Ad = NULL, adj; 134292bb6962SLisandro Dalcin IS ispart, isnumb, *is; 134392bb6962SLisandro Dalcin 134492bb6962SLisandro Dalcin PetscFunctionBegin; 13450700a824SBarry Smith PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 13464f572ea9SToby Isaac PetscAssertPointer(outis, 3); 134763a3b9bcSJacob Faibussowitsch PetscCheck(n >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "number of local blocks must be > 0, n = %" PetscInt_FMT, n); 134892bb6962SLisandro Dalcin 134992bb6962SLisandro Dalcin /* Get prefix, row distribution, and block size */ 13509566063dSJacob Faibussowitsch PetscCall(MatGetOptionsPrefix(A, &prefix)); 13519566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A, &rstart, &rend)); 13529566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(A, &bs)); 135363a3b9bcSJacob 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); 135465e19b50SBarry Smith 135592bb6962SLisandro Dalcin /* Get diagonal block from matrix if possible */ 13569566063dSJacob Faibussowitsch PetscCall(MatHasOperation(A, MATOP_GET_DIAGONAL_BLOCK, &hasop)); 135748a46eb9SPierre Jolivet if (hasop) PetscCall(MatGetDiagonalBlock(A, &Ad)); 135892bb6962SLisandro Dalcin if (Ad) { 13599566063dSJacob Faibussowitsch PetscCall(PetscObjectBaseTypeCompare((PetscObject)Ad, MATSEQBAIJ, &isbaij)); 13609566063dSJacob Faibussowitsch if (!isbaij) PetscCall(PetscObjectBaseTypeCompare((PetscObject)Ad, MATSEQSBAIJ, &isbaij)); 136192bb6962SLisandro Dalcin } 136292bb6962SLisandro Dalcin if (Ad && n > 1) { 1363ace3abfcSBarry Smith PetscBool match, done; 136492bb6962SLisandro Dalcin /* Try to setup a good matrix partitioning if available */ 13659566063dSJacob Faibussowitsch PetscCall(MatPartitioningCreate(PETSC_COMM_SELF, &mpart)); 13669566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)mpart, prefix)); 13679566063dSJacob Faibussowitsch PetscCall(MatPartitioningSetFromOptions(mpart)); 13689566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)mpart, MATPARTITIONINGCURRENT, &match)); 136948a46eb9SPierre Jolivet if (!match) PetscCall(PetscObjectTypeCompare((PetscObject)mpart, MATPARTITIONINGSQUARE, &match)); 137092bb6962SLisandro Dalcin if (!match) { /* assume a "good" partitioner is available */ 13711a83f524SJed Brown PetscInt na; 13721a83f524SJed Brown const PetscInt *ia, *ja; 13739566063dSJacob Faibussowitsch PetscCall(MatGetRowIJ(Ad, 0, PETSC_TRUE, isbaij, &na, &ia, &ja, &done)); 137492bb6962SLisandro Dalcin if (done) { 137592bb6962SLisandro Dalcin /* Build adjacency matrix by hand. Unfortunately a call to 137692bb6962SLisandro Dalcin MatConvert(Ad,MATMPIADJ,MAT_INITIAL_MATRIX,&adj) will 137792bb6962SLisandro Dalcin remove the block-aij structure and we cannot expect 137892bb6962SLisandro Dalcin MatPartitioning to split vertices as we need */ 13790a545947SLisandro Dalcin PetscInt i, j, len, nnz, cnt, *iia = NULL, *jja = NULL; 13801a83f524SJed Brown const PetscInt *row; 138192bb6962SLisandro Dalcin nnz = 0; 138292bb6962SLisandro Dalcin for (i = 0; i < na; i++) { /* count number of nonzeros */ 138392bb6962SLisandro Dalcin len = ia[i + 1] - ia[i]; 138492bb6962SLisandro Dalcin row = ja + ia[i]; 138592bb6962SLisandro Dalcin for (j = 0; j < len; j++) { 138692bb6962SLisandro Dalcin if (row[j] == i) { /* don't count diagonal */ 13879371c9d4SSatish Balay len--; 13889371c9d4SSatish Balay break; 138992bb6962SLisandro Dalcin } 139092bb6962SLisandro Dalcin } 139192bb6962SLisandro Dalcin nnz += len; 139292bb6962SLisandro Dalcin } 13939566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(na + 1, &iia)); 13949566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nnz, &jja)); 139592bb6962SLisandro Dalcin nnz = 0; 139692bb6962SLisandro Dalcin iia[0] = 0; 139792bb6962SLisandro Dalcin for (i = 0; i < na; i++) { /* fill adjacency */ 139892bb6962SLisandro Dalcin cnt = 0; 139992bb6962SLisandro Dalcin len = ia[i + 1] - ia[i]; 140092bb6962SLisandro Dalcin row = ja + ia[i]; 140192bb6962SLisandro Dalcin for (j = 0; j < len; j++) { 140292bb6962SLisandro Dalcin if (row[j] != i) { /* if not diagonal */ 140392bb6962SLisandro Dalcin jja[nnz + cnt++] = row[j]; 140492bb6962SLisandro Dalcin } 140592bb6962SLisandro Dalcin } 140692bb6962SLisandro Dalcin nnz += cnt; 140792bb6962SLisandro Dalcin iia[i + 1] = nnz; 140892bb6962SLisandro Dalcin } 140992bb6962SLisandro Dalcin /* Partitioning of the adjacency matrix */ 14109566063dSJacob Faibussowitsch PetscCall(MatCreateMPIAdj(PETSC_COMM_SELF, na, na, iia, jja, NULL, &adj)); 14119566063dSJacob Faibussowitsch PetscCall(MatPartitioningSetAdjacency(mpart, adj)); 14129566063dSJacob Faibussowitsch PetscCall(MatPartitioningSetNParts(mpart, n)); 14139566063dSJacob Faibussowitsch PetscCall(MatPartitioningApply(mpart, &ispart)); 14149566063dSJacob Faibussowitsch PetscCall(ISPartitioningToNumbering(ispart, &isnumb)); 14159566063dSJacob Faibussowitsch PetscCall(MatDestroy(&adj)); 141692bb6962SLisandro Dalcin foundpart = PETSC_TRUE; 141792bb6962SLisandro Dalcin } 14189566063dSJacob Faibussowitsch PetscCall(MatRestoreRowIJ(Ad, 0, PETSC_TRUE, isbaij, &na, &ia, &ja, &done)); 141992bb6962SLisandro Dalcin } 14209566063dSJacob Faibussowitsch PetscCall(MatPartitioningDestroy(&mpart)); 142192bb6962SLisandro Dalcin } 142292bb6962SLisandro Dalcin 14239566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &is)); 142492bb6962SLisandro Dalcin *outis = is; 142592bb6962SLisandro Dalcin 142692bb6962SLisandro Dalcin if (!foundpart) { 142792bb6962SLisandro Dalcin /* Partitioning by contiguous chunks of rows */ 142892bb6962SLisandro Dalcin 142992bb6962SLisandro Dalcin PetscInt mbs = (rend - rstart) / bs; 143092bb6962SLisandro Dalcin PetscInt start = rstart; 143192bb6962SLisandro Dalcin for (i = 0; i < n; i++) { 143292bb6962SLisandro Dalcin PetscInt count = (mbs / n + ((mbs % n) > i)) * bs; 14339566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, count, start, 1, &is[i])); 143492bb6962SLisandro Dalcin start += count; 143592bb6962SLisandro Dalcin } 143692bb6962SLisandro Dalcin 143792bb6962SLisandro Dalcin } else { 143892bb6962SLisandro Dalcin /* Partitioning by adjacency of diagonal block */ 143992bb6962SLisandro Dalcin 144092bb6962SLisandro Dalcin const PetscInt *numbering; 144192bb6962SLisandro Dalcin PetscInt *count, nidx, *indices, *newidx, start = 0; 144292bb6962SLisandro Dalcin /* Get node count in each partition */ 14439566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &count)); 14449566063dSJacob Faibussowitsch PetscCall(ISPartitioningCount(ispart, n, count)); 144592bb6962SLisandro Dalcin if (isbaij && bs > 1) { /* adjust for the block-aij case */ 144692bb6962SLisandro Dalcin for (i = 0; i < n; i++) count[i] *= bs; 144792bb6962SLisandro Dalcin } 144892bb6962SLisandro Dalcin /* Build indices from node numbering */ 14499566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(isnumb, &nidx)); 14509566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nidx, &indices)); 145192bb6962SLisandro Dalcin for (i = 0; i < nidx; i++) indices[i] = i; /* needs to be initialized */ 14529566063dSJacob Faibussowitsch PetscCall(ISGetIndices(isnumb, &numbering)); 14539566063dSJacob Faibussowitsch PetscCall(PetscSortIntWithPermutation(nidx, numbering, indices)); 14549566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(isnumb, &numbering)); 145592bb6962SLisandro Dalcin if (isbaij && bs > 1) { /* adjust for the block-aij case */ 14569566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nidx * bs, &newidx)); 14572fa5cd67SKarl Rupp for (i = 0; i < nidx; i++) { 14582fa5cd67SKarl Rupp for (j = 0; j < bs; j++) newidx[i * bs + j] = indices[i] * bs + j; 14592fa5cd67SKarl Rupp } 14609566063dSJacob Faibussowitsch PetscCall(PetscFree(indices)); 146192bb6962SLisandro Dalcin nidx *= bs; 146292bb6962SLisandro Dalcin indices = newidx; 146392bb6962SLisandro Dalcin } 146492bb6962SLisandro Dalcin /* Shift to get global indices */ 146592bb6962SLisandro Dalcin for (i = 0; i < nidx; i++) indices[i] += rstart; 146692bb6962SLisandro Dalcin 146792bb6962SLisandro Dalcin /* Build the index sets for each block */ 146892bb6962SLisandro Dalcin for (i = 0; i < n; i++) { 14699566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, count[i], &indices[start], PETSC_COPY_VALUES, &is[i])); 14709566063dSJacob Faibussowitsch PetscCall(ISSort(is[i])); 147192bb6962SLisandro Dalcin start += count[i]; 147292bb6962SLisandro Dalcin } 147392bb6962SLisandro Dalcin 14749566063dSJacob Faibussowitsch PetscCall(PetscFree(count)); 14759566063dSJacob Faibussowitsch PetscCall(PetscFree(indices)); 14769566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isnumb)); 14779566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ispart)); 147892bb6962SLisandro Dalcin } 14793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 148092bb6962SLisandro Dalcin } 148192bb6962SLisandro Dalcin 148292bb6962SLisandro Dalcin /*@C 148392bb6962SLisandro Dalcin PCASMDestroySubdomains - Destroys the index sets created with 1484f1580f4eSBarry Smith `PCASMCreateSubdomains()`. Should be called after setting subdomains with `PCASMSetLocalSubdomains()`. 148592bb6962SLisandro Dalcin 148692bb6962SLisandro Dalcin Collective 148792bb6962SLisandro Dalcin 148892bb6962SLisandro Dalcin Input Parameters: 148992bb6962SLisandro Dalcin + n - the number of index sets 14902b691e39SMatthew Knepley . is - the array of index sets 14912fe279fdSBarry Smith - is_local - the array of local index sets, can be `NULL` 149292bb6962SLisandro Dalcin 149392bb6962SLisandro Dalcin Level: advanced 149492bb6962SLisandro Dalcin 1495ce78bad3SBarry Smith Developer Note: 1496ce78bad3SBarry Smith The `IS` arguments should be a *[] 1497ce78bad3SBarry Smith 1498562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMCreateSubdomains()`, `PCASMSetLocalSubdomains()` 149992bb6962SLisandro Dalcin @*/ 1500ce78bad3SBarry Smith PetscErrorCode PCASMDestroySubdomains(PetscInt n, IS *is[], IS *is_local[]) 1501d71ae5a4SJacob Faibussowitsch { 150292bb6962SLisandro Dalcin PetscInt i; 15035fd66863SKarl Rupp 150492bb6962SLisandro Dalcin PetscFunctionBegin; 15053ba16761SJacob Faibussowitsch if (n <= 0) PetscFunctionReturn(PETSC_SUCCESS); 1506ce78bad3SBarry Smith if (*is) { 1507ce78bad3SBarry Smith PetscAssertPointer(*is, 2); 1508ce78bad3SBarry Smith for (i = 0; i < n; i++) PetscCall(ISDestroy(&(*is)[i])); 1509ce78bad3SBarry Smith PetscCall(PetscFree(*is)); 1510a35b7b57SDmitry Karpeev } 1511ce78bad3SBarry Smith if (is_local && *is_local) { 1512ce78bad3SBarry Smith PetscAssertPointer(*is_local, 3); 1513ce78bad3SBarry Smith for (i = 0; i < n; i++) PetscCall(ISDestroy(&(*is_local)[i])); 1514ce78bad3SBarry Smith PetscCall(PetscFree(*is_local)); 15152b691e39SMatthew Knepley } 15163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 151792bb6962SLisandro Dalcin } 151892bb6962SLisandro Dalcin 15196141accfSBarry Smith /*@C 15204b9ad928SBarry Smith PCASMCreateSubdomains2D - Creates the index sets for the overlapping Schwarz 1521f1580f4eSBarry Smith preconditioner, `PCASM`, for a two-dimensional problem on a regular grid. 15224b9ad928SBarry Smith 15234b9ad928SBarry Smith Not Collective 15244b9ad928SBarry Smith 15254b9ad928SBarry Smith Input Parameters: 15266b867d5aSJose E. Roman + m - the number of mesh points in the x direction 15276b867d5aSJose E. Roman . n - the number of mesh points in the y direction 15286b867d5aSJose E. Roman . M - the number of subdomains in the x direction 15296b867d5aSJose E. Roman . N - the number of subdomains in the y direction 15304b9ad928SBarry Smith . dof - degrees of freedom per node 15314b9ad928SBarry Smith - overlap - overlap in mesh lines 15324b9ad928SBarry Smith 15334b9ad928SBarry Smith Output Parameters: 15344b9ad928SBarry Smith + Nsub - the number of subdomains created 15353d03bb48SJed Brown . is - array of index sets defining overlapping (if overlap > 0) subdomains 15363d03bb48SJed Brown - is_local - array of index sets defining non-overlapping subdomains 15374b9ad928SBarry Smith 153820f4b53cSBarry Smith Level: advanced 153920f4b53cSBarry Smith 15404b9ad928SBarry Smith Note: 1541f1580f4eSBarry Smith Presently `PCAMSCreateSubdomains2d()` is valid only for sequential 15424b9ad928SBarry Smith preconditioners. More general related routines are 1543f1580f4eSBarry Smith `PCASMSetTotalSubdomains()` and `PCASMSetLocalSubdomains()`. 15444b9ad928SBarry Smith 1545562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMSetTotalSubdomains()`, `PCASMSetLocalSubdomains()`, `PCASMGetSubKSP()`, 1546db781477SPatrick Sanan `PCASMSetOverlap()` 15474b9ad928SBarry Smith @*/ 1548ce78bad3SBarry Smith PetscErrorCode PCASMCreateSubdomains2D(PetscInt m, PetscInt n, PetscInt M, PetscInt N, PetscInt dof, PetscInt overlap, PetscInt *Nsub, IS *is[], IS *is_local[]) 1549d71ae5a4SJacob Faibussowitsch { 15503d03bb48SJed Brown PetscInt i, j, height, width, ystart, xstart, yleft, yright, xleft, xright, loc_outer; 155113f74950SBarry Smith PetscInt nidx, *idx, loc, ii, jj, count; 15524b9ad928SBarry Smith 15534b9ad928SBarry Smith PetscFunctionBegin; 15547827d75bSBarry Smith PetscCheck(dof == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "dof must be 1"); 15554b9ad928SBarry Smith 15564b9ad928SBarry Smith *Nsub = N * M; 15579566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(*Nsub, is)); 15589566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(*Nsub, is_local)); 15594b9ad928SBarry Smith ystart = 0; 15603d03bb48SJed Brown loc_outer = 0; 15614b9ad928SBarry Smith for (i = 0; i < N; i++) { 15624b9ad928SBarry Smith height = n / N + ((n % N) > i); /* height of subdomain */ 15637827d75bSBarry Smith PetscCheck(height >= 2, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many N subdomains for mesh dimension n"); 15649371c9d4SSatish Balay yleft = ystart - overlap; 15659371c9d4SSatish Balay if (yleft < 0) yleft = 0; 15669371c9d4SSatish Balay yright = ystart + height + overlap; 15679371c9d4SSatish Balay if (yright > n) yright = n; 15684b9ad928SBarry Smith xstart = 0; 15694b9ad928SBarry Smith for (j = 0; j < M; j++) { 15704b9ad928SBarry Smith width = m / M + ((m % M) > j); /* width of subdomain */ 15717827d75bSBarry Smith PetscCheck(width >= 2, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many M subdomains for mesh dimension m"); 15729371c9d4SSatish Balay xleft = xstart - overlap; 15739371c9d4SSatish Balay if (xleft < 0) xleft = 0; 15749371c9d4SSatish Balay xright = xstart + width + overlap; 15759371c9d4SSatish Balay if (xright > m) xright = m; 15764b9ad928SBarry Smith nidx = (xright - xleft) * (yright - yleft); 15779566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nidx, &idx)); 15784b9ad928SBarry Smith loc = 0; 15794b9ad928SBarry Smith for (ii = yleft; ii < yright; ii++) { 15804b9ad928SBarry Smith count = m * ii + xleft; 15812fa5cd67SKarl Rupp for (jj = xleft; jj < xright; jj++) idx[loc++] = count++; 15824b9ad928SBarry Smith } 15839566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, nidx, idx, PETSC_COPY_VALUES, (*is) + loc_outer)); 15843d03bb48SJed Brown if (overlap == 0) { 15859566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)(*is)[loc_outer])); 15862fa5cd67SKarl Rupp 15873d03bb48SJed Brown (*is_local)[loc_outer] = (*is)[loc_outer]; 15883d03bb48SJed Brown } else { 15893d03bb48SJed Brown for (loc = 0, ii = ystart; ii < ystart + height; ii++) { 1590ad540459SPierre Jolivet for (jj = xstart; jj < xstart + width; jj++) idx[loc++] = m * ii + jj; 15913d03bb48SJed Brown } 15929566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, loc, idx, PETSC_COPY_VALUES, *is_local + loc_outer)); 15933d03bb48SJed Brown } 15949566063dSJacob Faibussowitsch PetscCall(PetscFree(idx)); 15954b9ad928SBarry Smith xstart += width; 15963d03bb48SJed Brown loc_outer++; 15974b9ad928SBarry Smith } 15984b9ad928SBarry Smith ystart += height; 15994b9ad928SBarry Smith } 16009566063dSJacob Faibussowitsch for (i = 0; i < *Nsub; i++) PetscCall(ISSort((*is)[i])); 16013ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 16024b9ad928SBarry Smith } 16034b9ad928SBarry Smith 16044b9ad928SBarry Smith /*@C 16054b9ad928SBarry Smith PCASMGetLocalSubdomains - Gets the local subdomains (for this processor 1606f1580f4eSBarry Smith only) for the additive Schwarz preconditioner, `PCASM`. 16074b9ad928SBarry Smith 1608ad4df100SBarry Smith Not Collective 16094b9ad928SBarry Smith 16104b9ad928SBarry Smith Input Parameter: 16114b9ad928SBarry Smith . pc - the preconditioner context 16124b9ad928SBarry Smith 16134b9ad928SBarry Smith Output Parameters: 1614f17a6784SPierre Jolivet + n - if requested, the number of subdomains for this processor (default value = 1) 1615f17a6784SPierre Jolivet . is - if requested, the index sets that define the subdomains for this processor 161620f4b53cSBarry Smith - is_local - if requested, the index sets that define the local part of the subdomains for this processor (can be `NULL`) 161720f4b53cSBarry Smith 161820f4b53cSBarry Smith Level: advanced 16194b9ad928SBarry Smith 1620f1580f4eSBarry Smith Note: 1621f1580f4eSBarry Smith The `IS` numbering is in the parallel, global numbering of the vector. 16224b9ad928SBarry Smith 1623562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMSetTotalSubdomains()`, `PCASMSetOverlap()`, `PCASMGetSubKSP()`, 1624db781477SPatrick Sanan `PCASMCreateSubdomains2D()`, `PCASMSetLocalSubdomains()`, `PCASMGetLocalSubmatrices()` 16254b9ad928SBarry Smith @*/ 1626d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMGetLocalSubdomains(PC pc, PetscInt *n, IS *is[], IS *is_local[]) 1627d71ae5a4SJacob Faibussowitsch { 16282a808120SBarry Smith PC_ASM *osm = (PC_ASM *)pc->data; 1629ace3abfcSBarry Smith PetscBool match; 16304b9ad928SBarry Smith 16314b9ad928SBarry Smith PetscFunctionBegin; 16320700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 16334f572ea9SToby Isaac if (n) PetscAssertPointer(n, 2); 16344f572ea9SToby Isaac if (is) PetscAssertPointer(is, 3); 16354f572ea9SToby Isaac if (is_local) PetscAssertPointer(is_local, 4); 16369566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCASM, &match)); 163728b400f6SJacob Faibussowitsch PetscCheck(match, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "PC is not a PCASM"); 16384b9ad928SBarry Smith if (n) *n = osm->n_local_true; 16394b9ad928SBarry Smith if (is) *is = osm->is; 16402b691e39SMatthew Knepley if (is_local) *is_local = osm->is_local; 16413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 16424b9ad928SBarry Smith } 16434b9ad928SBarry Smith 16444b9ad928SBarry Smith /*@C 16454b9ad928SBarry Smith PCASMGetLocalSubmatrices - Gets the local submatrices (for this processor 1646f1580f4eSBarry Smith only) for the additive Schwarz preconditioner, `PCASM`. 16474b9ad928SBarry Smith 1648ad4df100SBarry Smith Not Collective 16494b9ad928SBarry Smith 16504b9ad928SBarry Smith Input Parameter: 16514b9ad928SBarry Smith . pc - the preconditioner context 16524b9ad928SBarry Smith 16534b9ad928SBarry Smith Output Parameters: 1654f17a6784SPierre Jolivet + n - if requested, the number of matrices for this processor (default value = 1) 1655f17a6784SPierre Jolivet - mat - if requested, the matrices 16564b9ad928SBarry Smith 16574b9ad928SBarry Smith Level: advanced 16584b9ad928SBarry Smith 165995452b02SPatrick Sanan Notes: 1660f1580f4eSBarry Smith Call after `PCSetUp()` (or `KSPSetUp()`) but before `PCApply()` and before `PCSetUpOnBlocks()`) 1661cf739d55SBarry Smith 1662f1580f4eSBarry Smith Usually one would use `PCSetModifySubMatrices()` to change the submatrices in building the preconditioner. 1663cf739d55SBarry Smith 1664562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMSetTotalSubdomains()`, `PCASMSetOverlap()`, `PCASMGetSubKSP()`, 1665db781477SPatrick Sanan `PCASMCreateSubdomains2D()`, `PCASMSetLocalSubdomains()`, `PCASMGetLocalSubdomains()`, `PCSetModifySubMatrices()` 16664b9ad928SBarry Smith @*/ 1667d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMGetLocalSubmatrices(PC pc, PetscInt *n, Mat *mat[]) 1668d71ae5a4SJacob Faibussowitsch { 16694b9ad928SBarry Smith PC_ASM *osm; 1670ace3abfcSBarry Smith PetscBool match; 16714b9ad928SBarry Smith 16724b9ad928SBarry Smith PetscFunctionBegin; 16730700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 16744f572ea9SToby Isaac if (n) PetscAssertPointer(n, 2); 16754f572ea9SToby Isaac if (mat) PetscAssertPointer(mat, 3); 167628b400f6SJacob Faibussowitsch PetscCheck(pc->setupcalled, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must call after KSPSetUp() or PCSetUp()."); 16779566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCASM, &match)); 167892bb6962SLisandro Dalcin if (!match) { 167992bb6962SLisandro Dalcin if (n) *n = 0; 16800298fd71SBarry Smith if (mat) *mat = NULL; 168192bb6962SLisandro Dalcin } else { 16824b9ad928SBarry Smith osm = (PC_ASM *)pc->data; 16834b9ad928SBarry Smith if (n) *n = osm->n_local_true; 16844b9ad928SBarry Smith if (mat) *mat = osm->pmat; 168592bb6962SLisandro Dalcin } 16863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 16874b9ad928SBarry Smith } 1688d709ea83SDmitry Karpeev 1689d709ea83SDmitry Karpeev /*@ 1690f1580f4eSBarry Smith PCASMSetDMSubdomains - Indicates whether to use `DMCreateDomainDecomposition()` to define the subdomains, whenever possible. 1691f1ee410cSBarry Smith 1692d709ea83SDmitry Karpeev Logically Collective 1693d709ea83SDmitry Karpeev 1694d8d19677SJose E. Roman Input Parameters: 1695d709ea83SDmitry Karpeev + pc - the preconditioner 1696f1580f4eSBarry Smith - flg - boolean indicating whether to use subdomains defined by the `DM` 1697d709ea83SDmitry Karpeev 1698d709ea83SDmitry Karpeev Options Database Key: 169973ff1848SBarry Smith . -pc_asm_dm_subdomains <bool> - use subdomains defined by the `DM` with `DMCreateDomainDecomposition()` 1700d709ea83SDmitry Karpeev 1701d709ea83SDmitry Karpeev Level: intermediate 1702d709ea83SDmitry Karpeev 1703f1580f4eSBarry Smith Note: 1704f1580f4eSBarry Smith `PCASMSetTotalSubdomains()` and `PCASMSetOverlap()` take precedence over `PCASMSetDMSubdomains()`, 1705d709ea83SDmitry Karpeev so setting either of the first two effectively turns the latter off. 1706d709ea83SDmitry Karpeev 170773ff1848SBarry Smith Developer Note: 170873ff1848SBarry Smith This should be `PCASMSetUseDMSubdomains()`, similarly for the options database key 170973ff1848SBarry Smith 1710562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMGetDMSubdomains()`, `PCASMSetTotalSubdomains()`, `PCASMSetOverlap()` 1711db781477SPatrick Sanan `PCASMCreateSubdomains2D()`, `PCASMSetLocalSubdomains()`, `PCASMGetLocalSubdomains()` 1712d709ea83SDmitry Karpeev @*/ 1713d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMSetDMSubdomains(PC pc, PetscBool flg) 1714d71ae5a4SJacob Faibussowitsch { 1715d709ea83SDmitry Karpeev PC_ASM *osm = (PC_ASM *)pc->data; 1716d709ea83SDmitry Karpeev PetscBool match; 1717d709ea83SDmitry Karpeev 1718d709ea83SDmitry Karpeev PetscFunctionBegin; 1719d709ea83SDmitry Karpeev PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1720d709ea83SDmitry Karpeev PetscValidLogicalCollectiveBool(pc, flg, 2); 172128b400f6SJacob Faibussowitsch PetscCheck(!pc->setupcalled, ((PetscObject)pc)->comm, PETSC_ERR_ARG_WRONGSTATE, "Not for a setup PC."); 17229566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCASM, &match)); 1723ad540459SPierre Jolivet if (match) osm->dm_subdomains = flg; 17243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1725d709ea83SDmitry Karpeev } 1726d709ea83SDmitry Karpeev 1727d709ea83SDmitry Karpeev /*@ 1728f1580f4eSBarry Smith PCASMGetDMSubdomains - Returns flag indicating whether to use `DMCreateDomainDecomposition()` to define the subdomains, whenever possible. 1729f1580f4eSBarry Smith 1730d709ea83SDmitry Karpeev Not Collective 1731d709ea83SDmitry Karpeev 1732d709ea83SDmitry Karpeev Input Parameter: 1733d709ea83SDmitry Karpeev . pc - the preconditioner 1734d709ea83SDmitry Karpeev 1735d709ea83SDmitry Karpeev Output Parameter: 173620f4b53cSBarry Smith . flg - boolean indicating whether to use subdomains defined by the `DM` 1737d709ea83SDmitry Karpeev 1738d709ea83SDmitry Karpeev Level: intermediate 1739d709ea83SDmitry Karpeev 174073ff1848SBarry Smith Developer Note: 174173ff1848SBarry Smith This should be `PCASMSetUseDMSubdomains()` 174273ff1848SBarry Smith 1743562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMSetDMSubdomains()`, `PCASMSetTotalSubdomains()`, `PCASMSetOverlap()` 1744db781477SPatrick Sanan `PCASMCreateSubdomains2D()`, `PCASMSetLocalSubdomains()`, `PCASMGetLocalSubdomains()` 1745d709ea83SDmitry Karpeev @*/ 1746d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMGetDMSubdomains(PC pc, PetscBool *flg) 1747d71ae5a4SJacob Faibussowitsch { 1748d709ea83SDmitry Karpeev PC_ASM *osm = (PC_ASM *)pc->data; 1749d709ea83SDmitry Karpeev PetscBool match; 1750d709ea83SDmitry Karpeev 1751d709ea83SDmitry Karpeev PetscFunctionBegin; 1752d709ea83SDmitry Karpeev PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 17534f572ea9SToby Isaac PetscAssertPointer(flg, 2); 17549566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCASM, &match)); 175556ea09ceSStefano Zampini if (match) *flg = osm->dm_subdomains; 175656ea09ceSStefano Zampini else *flg = PETSC_FALSE; 17573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1758d709ea83SDmitry Karpeev } 175980ec0b7dSPatrick Sanan 17605d83a8b1SBarry Smith /*@ 1761f1580f4eSBarry Smith PCASMGetSubMatType - Gets the matrix type used for `PCASM` subsolves, as a string. 176280ec0b7dSPatrick Sanan 176380ec0b7dSPatrick Sanan Not Collective 176480ec0b7dSPatrick Sanan 176580ec0b7dSPatrick Sanan Input Parameter: 1766f1580f4eSBarry Smith . pc - the `PC` 176780ec0b7dSPatrick Sanan 176880ec0b7dSPatrick Sanan Output Parameter: 1769feefa0e1SJacob Faibussowitsch . sub_mat_type - name of matrix type 177080ec0b7dSPatrick Sanan 177180ec0b7dSPatrick Sanan Level: advanced 177280ec0b7dSPatrick Sanan 1773562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMSetSubMatType()`, `PCSetType()`, `VecSetType()`, `MatType`, `Mat` 177480ec0b7dSPatrick Sanan @*/ 1775d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMGetSubMatType(PC pc, MatType *sub_mat_type) 1776d71ae5a4SJacob Faibussowitsch { 177756ea09ceSStefano Zampini PetscFunctionBegin; 177856ea09ceSStefano Zampini PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1779cac4c232SBarry Smith PetscTryMethod(pc, "PCASMGetSubMatType_C", (PC, MatType *), (pc, sub_mat_type)); 17803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 178180ec0b7dSPatrick Sanan } 178280ec0b7dSPatrick Sanan 17835d83a8b1SBarry Smith /*@ 1784f1580f4eSBarry Smith PCASMSetSubMatType - Set the type of matrix used for `PCASM` subsolves 178580ec0b7dSPatrick Sanan 1786c3339decSBarry Smith Collective 178780ec0b7dSPatrick Sanan 178880ec0b7dSPatrick Sanan Input Parameters: 1789f1580f4eSBarry Smith + pc - the `PC` object 1790f1580f4eSBarry Smith - sub_mat_type - the `MatType` 179180ec0b7dSPatrick Sanan 179280ec0b7dSPatrick Sanan Options Database Key: 1793f1580f4eSBarry Smith . -pc_asm_sub_mat_type <sub_mat_type> - Sets the matrix type used for subsolves, for example, seqaijviennacl. 1794f1580f4eSBarry Smith If you specify a base name like aijviennacl, the corresponding sequential type is assumed. 179580ec0b7dSPatrick Sanan 1796f1580f4eSBarry Smith Note: 17972fe279fdSBarry Smith See `MatType` for available types 179880ec0b7dSPatrick Sanan 179980ec0b7dSPatrick Sanan Level: advanced 180080ec0b7dSPatrick Sanan 1801562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMGetSubMatType()`, `PCSetType()`, `VecSetType()`, `MatType`, `Mat` 180280ec0b7dSPatrick Sanan @*/ 1803d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMSetSubMatType(PC pc, MatType sub_mat_type) 1804d71ae5a4SJacob Faibussowitsch { 180556ea09ceSStefano Zampini PetscFunctionBegin; 180656ea09ceSStefano Zampini PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1807cac4c232SBarry Smith PetscTryMethod(pc, "PCASMSetSubMatType_C", (PC, MatType), (pc, sub_mat_type)); 18083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 180980ec0b7dSPatrick Sanan } 1810