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)); 5463a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " [%d] number of local blocks = %" PetscInt_FMT "\n", (int)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)); 62b4025f61SBarry Smith PetscCall(PetscViewerASCIIPrintf(sviewer, "[%d] local block number %" PetscInt_FMT ", size = %" PetscInt_FMT "\n", (int)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; 2071c2dc1cbSBarry Smith PetscCall(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; 5724b9ad928SBarry Smith /* 5734b9ad928SBarry Smith Support for limiting the restriction or interpolation to only local 5744b9ad928SBarry Smith subdomain values (leaving the other values 0). 5754b9ad928SBarry Smith 5764b9ad928SBarry Smith Note: these are reversed from the PCApply_ASM() because we are applying the 5774b9ad928SBarry Smith transpose of the three terms 5784b9ad928SBarry Smith */ 579d8b3b5e3Seaulisa 5804b9ad928SBarry Smith if (!(osm->type & PC_ASM_INTERPOLATE)) { 5814b9ad928SBarry Smith forward = SCATTER_FORWARD_LOCAL; 5824b9ad928SBarry Smith /* have to zero the work RHS since scatter may leave some slots empty */ 5839566063dSJacob Faibussowitsch PetscCall(VecSet(osm->lx, 0.0)); 5844b9ad928SBarry Smith } 5852fa5cd67SKarl Rupp if (!(osm->type & PC_ASM_RESTRICT)) reverse = SCATTER_REVERSE_LOCAL; 5864b9ad928SBarry Smith 587b0de9f23SBarry Smith /* zero the global and the local solutions */ 5889566063dSJacob Faibussowitsch PetscCall(VecSet(y, 0.0)); 5899566063dSJacob Faibussowitsch PetscCall(VecSet(osm->ly, 0.0)); 590d8b3b5e3Seaulisa 591b0de9f23SBarry Smith /* Copy the global RHS to local RHS including the ghost nodes */ 5929566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->restriction, x, osm->lx, INSERT_VALUES, forward)); 5939566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->restriction, x, osm->lx, INSERT_VALUES, forward)); 594d8b3b5e3Seaulisa 595b0de9f23SBarry Smith /* Restrict local RHS to the overlapping 0-block RHS */ 5969566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward)); 5979566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward)); 598d8b3b5e3Seaulisa 5994b9ad928SBarry Smith /* do the local solves */ 600d8b3b5e3Seaulisa for (i = 0; i < n_local_true; ++i) { 601b0de9f23SBarry Smith /* solve the overlapping i-block */ 6029566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PC_ApplyOnBlocks, osm->ksp[i], osm->x[i], osm->y[i], 0)); 6039566063dSJacob Faibussowitsch PetscCall(KSPSolveTranspose(osm->ksp[i], osm->x[i], osm->y[i])); 6049566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(osm->ksp[i], pc, osm->y[i])); 6059566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PC_ApplyOnBlocks, osm->ksp[i], osm->x[i], osm->y[i], 0)); 606d8b3b5e3Seaulisa 607a681a0f1SPierre Jolivet if (osm->lprolongation && osm->type != PC_ASM_RESTRICT) { /* interpolate the non-overlapping i-block solution to the local solution */ 6089566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward)); 6099566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward)); 610910cf402Sprj- } else { /* interpolate the overlapping i-block solution to the local solution */ 6119566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse)); 6129566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse)); 6134b9ad928SBarry Smith } 614d8b3b5e3Seaulisa 615d8b3b5e3Seaulisa if (i < n_local_true - 1) { 616b0de9f23SBarry Smith /* Restrict local RHS to the overlapping (i+1)-block RHS */ 6179566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->lrestriction[i + 1], osm->lx, osm->x[i + 1], INSERT_VALUES, forward)); 6189566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->lrestriction[i + 1], osm->lx, osm->x[i + 1], INSERT_VALUES, forward)); 6194b9ad928SBarry Smith } 6204b9ad928SBarry Smith } 621b0de9f23SBarry Smith /* Add the local solution to the global solution including the ghost nodes */ 6229566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->restriction, osm->ly, y, ADD_VALUES, reverse)); 6239566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->restriction, osm->ly, y, ADD_VALUES, reverse)); 6243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6254b9ad928SBarry Smith } 6264b9ad928SBarry Smith 627d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCReset_ASM(PC pc) 628d71ae5a4SJacob Faibussowitsch { 6294b9ad928SBarry Smith PC_ASM *osm = (PC_ASM *)pc->data; 63013f74950SBarry Smith PetscInt i; 6314b9ad928SBarry Smith 6324b9ad928SBarry Smith PetscFunctionBegin; 63392bb6962SLisandro Dalcin if (osm->ksp) { 63448a46eb9SPierre Jolivet for (i = 0; i < osm->n_local_true; i++) PetscCall(KSPReset(osm->ksp[i])); 63592bb6962SLisandro Dalcin } 636e09e08ccSBarry Smith if (osm->pmat) { 63748a46eb9SPierre Jolivet if (osm->n_local_true > 0) PetscCall(MatDestroySubMatrices(osm->n_local_true, &osm->pmat)); 63892bb6962SLisandro Dalcin } 6391dd8081eSeaulisa if (osm->lrestriction) { 6409566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&osm->restriction)); 6411dd8081eSeaulisa for (i = 0; i < osm->n_local_true; i++) { 6429566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&osm->lrestriction[i])); 6439566063dSJacob Faibussowitsch if (osm->lprolongation) PetscCall(VecScatterDestroy(&osm->lprolongation[i])); 6449566063dSJacob Faibussowitsch PetscCall(VecDestroy(&osm->x[i])); 6459566063dSJacob Faibussowitsch PetscCall(VecDestroy(&osm->y[i])); 6464b9ad928SBarry Smith } 6479566063dSJacob Faibussowitsch PetscCall(PetscFree(osm->lrestriction)); 6489566063dSJacob Faibussowitsch if (osm->lprolongation) PetscCall(PetscFree(osm->lprolongation)); 6499566063dSJacob Faibussowitsch PetscCall(PetscFree(osm->x)); 6509566063dSJacob Faibussowitsch PetscCall(PetscFree(osm->y)); 65192bb6962SLisandro Dalcin } 6529566063dSJacob Faibussowitsch PetscCall(PCASMDestroySubdomains(osm->n_local_true, osm->is, osm->is_local)); 6539566063dSJacob Faibussowitsch PetscCall(ISDestroy(&osm->lis)); 6549566063dSJacob Faibussowitsch PetscCall(VecDestroy(&osm->lx)); 6559566063dSJacob Faibussowitsch PetscCall(VecDestroy(&osm->ly)); 65648a46eb9SPierre Jolivet if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE) PetscCall(MatDestroyMatrices(osm->n_local_true, &osm->lmats)); 6572fa5cd67SKarl Rupp 6589566063dSJacob Faibussowitsch PetscCall(PetscFree(osm->sub_mat_type)); 65980ec0b7dSPatrick Sanan 6600a545947SLisandro Dalcin osm->is = NULL; 6610a545947SLisandro Dalcin osm->is_local = NULL; 6623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 663e91c6855SBarry Smith } 664e91c6855SBarry Smith 665d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_ASM(PC pc) 666d71ae5a4SJacob Faibussowitsch { 667e91c6855SBarry Smith PC_ASM *osm = (PC_ASM *)pc->data; 668e91c6855SBarry Smith PetscInt i; 669e91c6855SBarry Smith 670e91c6855SBarry Smith PetscFunctionBegin; 6719566063dSJacob Faibussowitsch PetscCall(PCReset_ASM(pc)); 672e91c6855SBarry Smith if (osm->ksp) { 67348a46eb9SPierre Jolivet for (i = 0; i < osm->n_local_true; i++) PetscCall(KSPDestroy(&osm->ksp[i])); 6749566063dSJacob Faibussowitsch PetscCall(PetscFree(osm->ksp)); 675e91c6855SBarry Smith } 6769566063dSJacob Faibussowitsch PetscCall(PetscFree(pc->data)); 67796322394SPierre Jolivet 6789566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetLocalSubdomains_C", NULL)); 6799566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetTotalSubdomains_C", NULL)); 6809566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetOverlap_C", NULL)); 6819566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetType_C", NULL)); 6829566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMGetType_C", NULL)); 6839566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetLocalType_C", NULL)); 6849566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMGetLocalType_C", NULL)); 6859566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetSortIndices_C", NULL)); 6869566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMGetSubKSP_C", NULL)); 6879566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMGetSubMatType_C", NULL)); 6889566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetSubMatType_C", NULL)); 6893ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6904b9ad928SBarry Smith } 6914b9ad928SBarry Smith 692d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetFromOptions_ASM(PC pc, PetscOptionItems *PetscOptionsObject) 693d71ae5a4SJacob Faibussowitsch { 6944b9ad928SBarry Smith PC_ASM *osm = (PC_ASM *)pc->data; 6959dcbbd2bSBarry Smith PetscInt blocks, ovl; 69657501b6eSBarry Smith PetscBool flg; 69792bb6962SLisandro Dalcin PCASMType asmtype; 69812cd4985SMatthew G. Knepley PCCompositeType loctype; 69980ec0b7dSPatrick Sanan char sub_mat_type[256]; 7004b9ad928SBarry Smith 7014b9ad928SBarry Smith PetscFunctionBegin; 702d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "Additive Schwarz options"); 7039566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_asm_dm_subdomains", "Use DMCreateDomainDecomposition() to define subdomains", "PCASMSetDMSubdomains", osm->dm_subdomains, &osm->dm_subdomains, &flg)); 7049566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_asm_blocks", "Number of subdomains", "PCASMSetTotalSubdomains", osm->n, &blocks, &flg)); 70565db9045SDmitry Karpeev if (flg) { 7069566063dSJacob Faibussowitsch PetscCall(PCASMSetTotalSubdomains(pc, blocks, NULL, NULL)); 707d709ea83SDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 70865db9045SDmitry Karpeev } 7099566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_asm_local_blocks", "Number of local subdomains", "PCASMSetLocalSubdomains", osm->n_local_true, &blocks, &flg)); 710342c94f9SMatthew G. Knepley if (flg) { 7119566063dSJacob Faibussowitsch PetscCall(PCASMSetLocalSubdomains(pc, blocks, NULL, NULL)); 712342c94f9SMatthew G. Knepley osm->dm_subdomains = PETSC_FALSE; 713342c94f9SMatthew G. Knepley } 7149566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_asm_overlap", "Number of grid points overlap", "PCASMSetOverlap", osm->overlap, &ovl, &flg)); 71565db9045SDmitry Karpeev if (flg) { 7169566063dSJacob Faibussowitsch PetscCall(PCASMSetOverlap(pc, ovl)); 717d709ea83SDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 71865db9045SDmitry Karpeev } 71990d69ab7SBarry Smith flg = PETSC_FALSE; 7209566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-pc_asm_type", "Type of restriction/extension", "PCASMSetType", PCASMTypes, (PetscEnum)osm->type, (PetscEnum *)&asmtype, &flg)); 7219566063dSJacob Faibussowitsch if (flg) PetscCall(PCASMSetType(pc, asmtype)); 72212cd4985SMatthew G. Knepley flg = PETSC_FALSE; 7239566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-pc_asm_local_type", "Type of local solver composition", "PCASMSetLocalType", PCCompositeTypes, (PetscEnum)osm->loctype, (PetscEnum *)&loctype, &flg)); 7249566063dSJacob Faibussowitsch if (flg) PetscCall(PCASMSetLocalType(pc, loctype)); 7259566063dSJacob Faibussowitsch PetscCall(PetscOptionsFList("-pc_asm_sub_mat_type", "Subsolve Matrix Type", "PCASMSetSubMatType", MatList, NULL, sub_mat_type, 256, &flg)); 7261baa6e33SBarry Smith if (flg) PetscCall(PCASMSetSubMatType(pc, sub_mat_type)); 727d0609cedSBarry Smith PetscOptionsHeadEnd(); 7283ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7294b9ad928SBarry Smith } 7304b9ad928SBarry Smith 731d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCASMSetLocalSubdomains_ASM(PC pc, PetscInt n, IS is[], IS is_local[]) 732d71ae5a4SJacob Faibussowitsch { 7334b9ad928SBarry Smith PC_ASM *osm = (PC_ASM *)pc->data; 73492bb6962SLisandro Dalcin PetscInt i; 7354b9ad928SBarry Smith 7364b9ad928SBarry Smith PetscFunctionBegin; 73763a3b9bcSJacob Faibussowitsch PetscCheck(n >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Each process must have 1 or more blocks, n = %" PetscInt_FMT, n); 7387827d75bSBarry Smith PetscCheck(!pc->setupcalled || (n == osm->n_local_true && !is), PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "PCASMSetLocalSubdomains() should be called before calling PCSetUp()."); 739e7e72b3dSBarry Smith 7404b9ad928SBarry Smith if (!pc->setupcalled) { 74192bb6962SLisandro Dalcin if (is) { 7429566063dSJacob Faibussowitsch for (i = 0; i < n; i++) PetscCall(PetscObjectReference((PetscObject)is[i])); 74392bb6962SLisandro Dalcin } 744832fc9a5SMatthew Knepley if (is_local) { 7459566063dSJacob Faibussowitsch for (i = 0; i < n; i++) PetscCall(PetscObjectReference((PetscObject)is_local[i])); 746832fc9a5SMatthew Knepley } 7479566063dSJacob Faibussowitsch PetscCall(PCASMDestroySubdomains(osm->n_local_true, osm->is, osm->is_local)); 7482fa5cd67SKarl Rupp 74907517c86SMark Adams if (osm->ksp && osm->n_local_true != n) { 75007517c86SMark Adams for (i = 0; i < osm->n_local_true; i++) PetscCall(KSPDestroy(&osm->ksp[i])); 75107517c86SMark Adams PetscCall(PetscFree(osm->ksp)); 75207517c86SMark Adams } 75307517c86SMark Adams 7544b9ad928SBarry Smith osm->n_local_true = n; 7550a545947SLisandro Dalcin osm->is = NULL; 7560a545947SLisandro Dalcin osm->is_local = NULL; 75792bb6962SLisandro Dalcin if (is) { 7589566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &osm->is)); 7592fa5cd67SKarl Rupp for (i = 0; i < n; i++) osm->is[i] = is[i]; 7603d03bb48SJed Brown /* Flag indicating that the user has set overlapping subdomains so PCASM should not increase their size. */ 7613d03bb48SJed Brown osm->overlap = -1; 76292bb6962SLisandro Dalcin } 7632b691e39SMatthew Knepley if (is_local) { 7649566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &osm->is_local)); 7652fa5cd67SKarl Rupp for (i = 0; i < n; i++) osm->is_local[i] = is_local[i]; 766a35b7b57SDmitry Karpeev if (!is) { 7679566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(osm->n_local_true, &osm->is)); 768a35b7b57SDmitry Karpeev for (i = 0; i < osm->n_local_true; i++) { 769a35b7b57SDmitry Karpeev if (osm->overlap > 0) { /* With positive overlap, osm->is[i] will be modified */ 7709566063dSJacob Faibussowitsch PetscCall(ISDuplicate(osm->is_local[i], &osm->is[i])); 7719566063dSJacob Faibussowitsch PetscCall(ISCopy(osm->is_local[i], osm->is[i])); 772a35b7b57SDmitry Karpeev } else { 7739566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)osm->is_local[i])); 774a35b7b57SDmitry Karpeev osm->is[i] = osm->is_local[i]; 775a35b7b57SDmitry Karpeev } 776a35b7b57SDmitry Karpeev } 777a35b7b57SDmitry Karpeev } 7782b691e39SMatthew Knepley } 7794b9ad928SBarry Smith } 7803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7814b9ad928SBarry Smith } 7824b9ad928SBarry Smith 783d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCASMSetTotalSubdomains_ASM(PC pc, PetscInt N, IS *is, IS *is_local) 784d71ae5a4SJacob Faibussowitsch { 7854b9ad928SBarry Smith PC_ASM *osm = (PC_ASM *)pc->data; 78613f74950SBarry Smith PetscMPIInt rank, size; 78778904715SLisandro Dalcin PetscInt n; 7884b9ad928SBarry Smith 7894b9ad928SBarry Smith PetscFunctionBegin; 79063a3b9bcSJacob Faibussowitsch PetscCheck(N >= 1, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Number of total blocks must be > 0, N = %" PetscInt_FMT, N); 79100045ab3SPierre Jolivet PetscCheck(!is && !is_local, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Use PCASMSetLocalSubdomains() to set specific index sets, they cannot be set globally yet."); 7924b9ad928SBarry Smith 7934b9ad928SBarry Smith /* 794880770e9SJed Brown Split the subdomains equally among all processors 7954b9ad928SBarry Smith */ 7969566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank)); 7979566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size)); 7984b9ad928SBarry Smith n = N / size + ((N % size) > rank); 79963a3b9bcSJacob Faibussowitsch PetscCheck(n, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Process %d must have at least one block: total processors %d total blocks %" PetscInt_FMT, (int)rank, (int)size, N); 8007827d75bSBarry Smith PetscCheck(!pc->setupcalled || n == osm->n_local_true, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "PCASMSetTotalSubdomains() should be called before PCSetUp()."); 8014b9ad928SBarry Smith if (!pc->setupcalled) { 8029566063dSJacob Faibussowitsch PetscCall(PCASMDestroySubdomains(osm->n_local_true, osm->is, osm->is_local)); 8032fa5cd67SKarl Rupp 8044b9ad928SBarry Smith osm->n_local_true = n; 8050a545947SLisandro Dalcin osm->is = NULL; 8060a545947SLisandro Dalcin osm->is_local = NULL; 8074b9ad928SBarry Smith } 8083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8094b9ad928SBarry Smith } 8104b9ad928SBarry Smith 811d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCASMSetOverlap_ASM(PC pc, PetscInt ovl) 812d71ae5a4SJacob Faibussowitsch { 81392bb6962SLisandro Dalcin PC_ASM *osm = (PC_ASM *)pc->data; 8144b9ad928SBarry Smith 8154b9ad928SBarry Smith PetscFunctionBegin; 8167827d75bSBarry Smith PetscCheck(ovl >= 0, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Negative overlap value requested"); 8177827d75bSBarry Smith PetscCheck(!pc->setupcalled || ovl == osm->overlap, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "PCASMSetOverlap() should be called before PCSetUp()."); 8182fa5cd67SKarl Rupp if (!pc->setupcalled) osm->overlap = ovl; 8193ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8204b9ad928SBarry Smith } 8214b9ad928SBarry Smith 822d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCASMSetType_ASM(PC pc, PCASMType type) 823d71ae5a4SJacob Faibussowitsch { 82492bb6962SLisandro Dalcin PC_ASM *osm = (PC_ASM *)pc->data; 8254b9ad928SBarry Smith 8264b9ad928SBarry Smith PetscFunctionBegin; 8274b9ad928SBarry Smith osm->type = type; 828bf108f30SBarry Smith osm->type_set = PETSC_TRUE; 8293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8304b9ad928SBarry Smith } 8314b9ad928SBarry Smith 832d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCASMGetType_ASM(PC pc, PCASMType *type) 833d71ae5a4SJacob Faibussowitsch { 834c60c7ad4SBarry Smith PC_ASM *osm = (PC_ASM *)pc->data; 835c60c7ad4SBarry Smith 836c60c7ad4SBarry Smith PetscFunctionBegin; 837c60c7ad4SBarry Smith *type = osm->type; 8383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 839c60c7ad4SBarry Smith } 840c60c7ad4SBarry Smith 841d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCASMSetLocalType_ASM(PC pc, PCCompositeType type) 842d71ae5a4SJacob Faibussowitsch { 84312cd4985SMatthew G. Knepley PC_ASM *osm = (PC_ASM *)pc->data; 84412cd4985SMatthew G. Knepley 84512cd4985SMatthew G. Knepley PetscFunctionBegin; 8467827d75bSBarry 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"); 84712cd4985SMatthew G. Knepley osm->loctype = type; 8483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 84912cd4985SMatthew G. Knepley } 85012cd4985SMatthew G. Knepley 851d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCASMGetLocalType_ASM(PC pc, PCCompositeType *type) 852d71ae5a4SJacob Faibussowitsch { 85312cd4985SMatthew G. Knepley PC_ASM *osm = (PC_ASM *)pc->data; 85412cd4985SMatthew G. Knepley 85512cd4985SMatthew G. Knepley PetscFunctionBegin; 85612cd4985SMatthew G. Knepley *type = osm->loctype; 8573ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 85812cd4985SMatthew G. Knepley } 85912cd4985SMatthew G. Knepley 860d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCASMSetSortIndices_ASM(PC pc, PetscBool doSort) 861d71ae5a4SJacob Faibussowitsch { 8626ed231c7SMatthew Knepley PC_ASM *osm = (PC_ASM *)pc->data; 8636ed231c7SMatthew Knepley 8646ed231c7SMatthew Knepley PetscFunctionBegin; 8656ed231c7SMatthew Knepley osm->sort_indices = doSort; 8663ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8676ed231c7SMatthew Knepley } 8686ed231c7SMatthew Knepley 869d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCASMGetSubKSP_ASM(PC pc, PetscInt *n_local, PetscInt *first_local, KSP **ksp) 870d71ae5a4SJacob Faibussowitsch { 87192bb6962SLisandro Dalcin PC_ASM *osm = (PC_ASM *)pc->data; 8724b9ad928SBarry Smith 8734b9ad928SBarry Smith PetscFunctionBegin; 8747827d75bSBarry 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"); 8754b9ad928SBarry Smith 8762fa5cd67SKarl Rupp if (n_local) *n_local = osm->n_local_true; 87792bb6962SLisandro Dalcin if (first_local) { 8789566063dSJacob Faibussowitsch PetscCallMPI(MPI_Scan(&osm->n_local_true, first_local, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)pc))); 87992bb6962SLisandro Dalcin *first_local -= osm->n_local_true; 88092bb6962SLisandro Dalcin } 881ed155784SPierre Jolivet if (ksp) *ksp = osm->ksp; 8823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8834b9ad928SBarry Smith } 8844b9ad928SBarry Smith 885d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCASMGetSubMatType_ASM(PC pc, MatType *sub_mat_type) 886d71ae5a4SJacob Faibussowitsch { 88780ec0b7dSPatrick Sanan PC_ASM *osm = (PC_ASM *)pc->data; 88880ec0b7dSPatrick Sanan 88980ec0b7dSPatrick Sanan PetscFunctionBegin; 89080ec0b7dSPatrick Sanan PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 8914f572ea9SToby Isaac PetscAssertPointer(sub_mat_type, 2); 89280ec0b7dSPatrick Sanan *sub_mat_type = osm->sub_mat_type; 8933ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 89480ec0b7dSPatrick Sanan } 89580ec0b7dSPatrick Sanan 896d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCASMSetSubMatType_ASM(PC pc, MatType sub_mat_type) 897d71ae5a4SJacob Faibussowitsch { 89880ec0b7dSPatrick Sanan PC_ASM *osm = (PC_ASM *)pc->data; 89980ec0b7dSPatrick Sanan 90080ec0b7dSPatrick Sanan PetscFunctionBegin; 90180ec0b7dSPatrick Sanan PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 9029566063dSJacob Faibussowitsch PetscCall(PetscFree(osm->sub_mat_type)); 9039566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(sub_mat_type, (char **)&osm->sub_mat_type)); 9043ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 90580ec0b7dSPatrick Sanan } 90680ec0b7dSPatrick Sanan 9074b9ad928SBarry Smith /*@C 908f1580f4eSBarry Smith PCASMSetLocalSubdomains - Sets the local subdomains (for this processor only) for the additive Schwarz preconditioner `PCASM`. 9094b9ad928SBarry Smith 910c3339decSBarry Smith Collective 9114b9ad928SBarry Smith 9124b9ad928SBarry Smith Input Parameters: 9134b9ad928SBarry Smith + pc - the preconditioner context 9144b9ad928SBarry Smith . n - the number of subdomains for this processor (default value = 1) 915*a3b724e8SBarry Smith . is - the index set that defines the subdomains for this processor (or `NULL` for PETSc to determine subdomains) 916f1ee410cSBarry 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 91720f4b53cSBarry Smith (or `NULL` to not provide these) 9184b9ad928SBarry Smith 919342c94f9SMatthew G. Knepley Options Database Key: 92020f4b53cSBarry Smith . -pc_asm_local_blocks <blks> - Sets number of local blocks 92120f4b53cSBarry Smith 92220f4b53cSBarry Smith Level: advanced 923342c94f9SMatthew G. Knepley 9244b9ad928SBarry Smith Notes: 925*a3b724e8SBarry Smith The `IS` numbering is in the parallel, global numbering of the vector for both `is` and `is_local` 9264b9ad928SBarry Smith 927f1580f4eSBarry Smith By default the `PCASM` preconditioner uses 1 block per processor. 9284b9ad928SBarry Smith 929f1580f4eSBarry Smith Use `PCASMSetTotalSubdomains()` to set the subdomains for all processors. 9304b9ad928SBarry Smith 931*a3b724e8SBarry Smith If `is_local` is provided and `PCASMType` is `PC_ASM_RESTRICT` then the solution only over the `is_local` region is interpolated 932*a3b724e8SBarry Smith back to form the global solution (this is the standard restricted additive Schwarz method, RASM) 933f1ee410cSBarry Smith 934*a3b724e8SBarry Smith If `is_local` is provided and `PCASMType` is `PC_ASM_INTERPOLATE` or `PC_ASM_NONE` then an error is generated since there is 935f1ee410cSBarry Smith no code to handle that case. 936f1ee410cSBarry Smith 937562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMSetTotalSubdomains()`, `PCASMSetOverlap()`, `PCASMGetSubKSP()`, 938f1580f4eSBarry Smith `PCASMCreateSubdomains2D()`, `PCASMGetLocalSubdomains()`, `PCASMType`, `PCASMSetType()`, `PCGASM` 9394b9ad928SBarry Smith @*/ 940d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMSetLocalSubdomains(PC pc, PetscInt n, IS is[], IS is_local[]) 941d71ae5a4SJacob Faibussowitsch { 9424b9ad928SBarry Smith PetscFunctionBegin; 9430700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 944cac4c232SBarry Smith PetscTryMethod(pc, "PCASMSetLocalSubdomains_C", (PC, PetscInt, IS[], IS[]), (pc, n, is, is_local)); 9453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9464b9ad928SBarry Smith } 9474b9ad928SBarry Smith 9484b9ad928SBarry Smith /*@C 949feb221f8SDmitry Karpeev PCASMSetTotalSubdomains - Sets the subdomains for all processors for the 950f1580f4eSBarry Smith additive Schwarz preconditioner, `PCASM`. 9514b9ad928SBarry Smith 952c3339decSBarry Smith Collective, all MPI ranks must pass in the same array of `IS` 9534b9ad928SBarry Smith 9544b9ad928SBarry Smith Input Parameters: 9554b9ad928SBarry Smith + pc - the preconditioner context 956feb221f8SDmitry Karpeev . N - the number of subdomains for all processors 957*a3b724e8SBarry Smith . is - the index sets that define the subdomains for all processors (or `NULL` to ask PETSc to determine the subdomains) 958*a3b724e8SBarry Smith - is_local - the index sets that define the local part of the subdomains for this processor (or `NULL` to not provide this information) 9594b9ad928SBarry Smith 9604b9ad928SBarry Smith Options Database Key: 9614b9ad928SBarry Smith . -pc_asm_blocks <blks> - Sets total blocks 9624b9ad928SBarry Smith 96320f4b53cSBarry Smith Level: advanced 96420f4b53cSBarry Smith 9654b9ad928SBarry Smith Notes: 966*a3b724e8SBarry Smith Currently you cannot use this to set the actual subdomains with the argument `is` or `is_local`. 9674b9ad928SBarry Smith 968f1580f4eSBarry Smith By default the `PCASM` preconditioner uses 1 block per processor. 9694b9ad928SBarry Smith 9704b9ad928SBarry Smith These index sets cannot be destroyed until after completion of the 971f1580f4eSBarry Smith linear solves for which the `PCASM` preconditioner is being used. 9724b9ad928SBarry Smith 973f1580f4eSBarry Smith Use `PCASMSetLocalSubdomains()` to set local subdomains. 9744b9ad928SBarry Smith 975f1580f4eSBarry Smith The `IS` numbering is in the parallel, global numbering of the vector for both is and is_local 9761093a601SBarry Smith 977562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMSetLocalSubdomains()`, `PCASMSetOverlap()`, `PCASMGetSubKSP()`, 978f1580f4eSBarry Smith `PCASMCreateSubdomains2D()`, `PCGASM` 9794b9ad928SBarry Smith @*/ 980d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMSetTotalSubdomains(PC pc, PetscInt N, IS is[], IS is_local[]) 981d71ae5a4SJacob Faibussowitsch { 9824b9ad928SBarry Smith PetscFunctionBegin; 9830700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 984cac4c232SBarry Smith PetscTryMethod(pc, "PCASMSetTotalSubdomains_C", (PC, PetscInt, IS[], IS[]), (pc, N, is, is_local)); 9853ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9864b9ad928SBarry Smith } 9874b9ad928SBarry Smith 9884b9ad928SBarry Smith /*@ 9894b9ad928SBarry Smith PCASMSetOverlap - Sets the overlap between a pair of subdomains for the 990f1580f4eSBarry Smith additive Schwarz preconditioner, `PCASM`. 9914b9ad928SBarry Smith 992c3339decSBarry Smith Logically Collective 9934b9ad928SBarry Smith 9944b9ad928SBarry Smith Input Parameters: 9954b9ad928SBarry Smith + pc - the preconditioner context 9964b9ad928SBarry Smith - ovl - the amount of overlap between subdomains (ovl >= 0, default value = 1) 9974b9ad928SBarry Smith 9984b9ad928SBarry Smith Options Database Key: 9994b9ad928SBarry Smith . -pc_asm_overlap <ovl> - Sets overlap 10004b9ad928SBarry Smith 100120f4b53cSBarry Smith Level: intermediate 100220f4b53cSBarry Smith 10034b9ad928SBarry Smith Notes: 1004f1580f4eSBarry Smith By default the `PCASM` preconditioner uses 1 block per processor. To use 1005f1580f4eSBarry Smith multiple blocks per perocessor, see `PCASMSetTotalSubdomains()` and 1006f1580f4eSBarry Smith `PCASMSetLocalSubdomains()` (and the option -pc_asm_blocks <blks>). 10074b9ad928SBarry Smith 10084b9ad928SBarry Smith The overlap defaults to 1, so if one desires that no additional 10094b9ad928SBarry Smith overlap be computed beyond what may have been set with a call to 1010f1580f4eSBarry Smith `PCASMSetTotalSubdomains()` or `PCASMSetLocalSubdomains()`, then ovl 10114b9ad928SBarry Smith must be set to be 0. In particular, if one does not explicitly set 10124b9ad928SBarry Smith the subdomains an application code, then all overlap would be computed 1013f1580f4eSBarry Smith internally by PETSc, and using an overlap of 0 would result in an `PCASM` 10144b9ad928SBarry Smith variant that is equivalent to the block Jacobi preconditioner. 10154b9ad928SBarry Smith 1016f1ee410cSBarry Smith The default algorithm used by PETSc to increase overlap is fast, but not scalable, 1017f1ee410cSBarry Smith use the option -mat_increase_overlap_scalable when the problem and number of processes is large. 1018f1ee410cSBarry Smith 10192fe279fdSBarry Smith One can define initial index sets with any overlap via 1020f1580f4eSBarry Smith `PCASMSetLocalSubdomains()`; the routine 1021f1580f4eSBarry Smith `PCASMSetOverlap()` merely allows PETSc to extend that overlap further 10224b9ad928SBarry Smith if desired. 10234b9ad928SBarry Smith 1024562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMSetTotalSubdomains()`, `PCASMSetLocalSubdomains()`, `PCASMGetSubKSP()`, 1025f1580f4eSBarry Smith `PCASMCreateSubdomains2D()`, `PCASMGetLocalSubdomains()`, `MatIncreaseOverlap()`, `PCGASM` 10264b9ad928SBarry Smith @*/ 1027d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMSetOverlap(PC pc, PetscInt ovl) 1028d71ae5a4SJacob Faibussowitsch { 10294b9ad928SBarry Smith PetscFunctionBegin; 10300700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1031c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc, ovl, 2); 1032cac4c232SBarry Smith PetscTryMethod(pc, "PCASMSetOverlap_C", (PC, PetscInt), (pc, ovl)); 10333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 10344b9ad928SBarry Smith } 10354b9ad928SBarry Smith 10364b9ad928SBarry Smith /*@ 10374b9ad928SBarry Smith PCASMSetType - Sets the type of restriction and interpolation used 1038f1580f4eSBarry Smith for local problems in the additive Schwarz method, `PCASM`. 10394b9ad928SBarry Smith 1040c3339decSBarry Smith Logically Collective 10414b9ad928SBarry Smith 10424b9ad928SBarry Smith Input Parameters: 10434b9ad928SBarry Smith + pc - the preconditioner context 1044f1580f4eSBarry Smith - type - variant of `PCASM`, one of 10454b9ad928SBarry Smith .vb 10464b9ad928SBarry Smith PC_ASM_BASIC - full interpolation and restriction 104782b5ce2aSStefano Zampini PC_ASM_RESTRICT - full restriction, local processor interpolation (default) 10484b9ad928SBarry Smith PC_ASM_INTERPOLATE - full interpolation, local processor restriction 10494b9ad928SBarry Smith PC_ASM_NONE - local processor restriction and interpolation 10504b9ad928SBarry Smith .ve 10514b9ad928SBarry Smith 10524b9ad928SBarry Smith Options Database Key: 1053f1580f4eSBarry Smith . -pc_asm_type [basic,restrict,interpolate,none] - Sets `PCASMType` 10544b9ad928SBarry Smith 105520f4b53cSBarry Smith Level: intermediate 105620f4b53cSBarry Smith 1057f1580f4eSBarry Smith Note: 1058f1580f4eSBarry Smith if the is_local arguments are passed to `PCASMSetLocalSubdomains()` then they are used when `PC_ASM_RESTRICT` has been selected 1059f1ee410cSBarry Smith to limit the local processor interpolation 1060f1ee410cSBarry Smith 1061562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMSetTotalSubdomains()`, `PCASMGetSubKSP()`, 1062f1580f4eSBarry Smith `PCASMCreateSubdomains2D()`, `PCASMType`, `PCASMSetLocalType()`, `PCASMGetLocalType()`, `PCGASM` 10634b9ad928SBarry Smith @*/ 1064d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMSetType(PC pc, PCASMType type) 1065d71ae5a4SJacob Faibussowitsch { 10664b9ad928SBarry Smith PetscFunctionBegin; 10670700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1068c5eb9154SBarry Smith PetscValidLogicalCollectiveEnum(pc, type, 2); 1069cac4c232SBarry Smith PetscTryMethod(pc, "PCASMSetType_C", (PC, PCASMType), (pc, type)); 10703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 10714b9ad928SBarry Smith } 10724b9ad928SBarry Smith 1073c60c7ad4SBarry Smith /*@ 1074c60c7ad4SBarry Smith PCASMGetType - Gets the type of restriction and interpolation used 1075f1580f4eSBarry Smith for local problems in the additive Schwarz method, `PCASM`. 1076c60c7ad4SBarry Smith 1077c3339decSBarry Smith Logically Collective 1078c60c7ad4SBarry Smith 1079c60c7ad4SBarry Smith Input Parameter: 1080c60c7ad4SBarry Smith . pc - the preconditioner context 1081c60c7ad4SBarry Smith 1082c60c7ad4SBarry Smith Output Parameter: 1083f1580f4eSBarry Smith . type - variant of `PCASM`, one of 1084c60c7ad4SBarry Smith .vb 1085c60c7ad4SBarry Smith PC_ASM_BASIC - full interpolation and restriction 1086c60c7ad4SBarry Smith PC_ASM_RESTRICT - full restriction, local processor interpolation 1087c60c7ad4SBarry Smith PC_ASM_INTERPOLATE - full interpolation, local processor restriction 1088c60c7ad4SBarry Smith PC_ASM_NONE - local processor restriction and interpolation 1089c60c7ad4SBarry Smith .ve 1090c60c7ad4SBarry Smith 1091c60c7ad4SBarry Smith Options Database Key: 1092f1580f4eSBarry Smith . -pc_asm_type [basic,restrict,interpolate,none] - Sets `PCASM` type 1093c60c7ad4SBarry Smith 1094c60c7ad4SBarry Smith Level: intermediate 1095c60c7ad4SBarry Smith 1096562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMSetTotalSubdomains()`, `PCASMGetSubKSP()`, `PCGASM`, 1097db781477SPatrick Sanan `PCASMCreateSubdomains2D()`, `PCASMType`, `PCASMSetType()`, `PCASMSetLocalType()`, `PCASMGetLocalType()` 1098c60c7ad4SBarry Smith @*/ 1099d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMGetType(PC pc, PCASMType *type) 1100d71ae5a4SJacob Faibussowitsch { 1101c60c7ad4SBarry Smith PetscFunctionBegin; 1102c60c7ad4SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1103cac4c232SBarry Smith PetscUseMethod(pc, "PCASMGetType_C", (PC, PCASMType *), (pc, type)); 11043ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1105c60c7ad4SBarry Smith } 1106c60c7ad4SBarry Smith 110712cd4985SMatthew G. Knepley /*@ 1108f1580f4eSBarry Smith PCASMSetLocalType - Sets the type of composition used for local problems in the additive Schwarz method, `PCASM`. 110912cd4985SMatthew G. Knepley 1110c3339decSBarry Smith Logically Collective 111112cd4985SMatthew G. Knepley 111212cd4985SMatthew G. Knepley Input Parameters: 111312cd4985SMatthew G. Knepley + pc - the preconditioner context 111412cd4985SMatthew G. Knepley - type - type of composition, one of 111512cd4985SMatthew G. Knepley .vb 111612cd4985SMatthew G. Knepley PC_COMPOSITE_ADDITIVE - local additive combination 111712cd4985SMatthew G. Knepley PC_COMPOSITE_MULTIPLICATIVE - local multiplicative combination 111812cd4985SMatthew G. Knepley .ve 111912cd4985SMatthew G. Knepley 112012cd4985SMatthew G. Knepley Options Database Key: 112112cd4985SMatthew G. Knepley . -pc_asm_local_type [additive,multiplicative] - Sets local solver composition type 112212cd4985SMatthew G. Knepley 112312cd4985SMatthew G. Knepley Level: intermediate 112412cd4985SMatthew G. Knepley 1125562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMSetType()`, `PCASMGetType()`, `PCASMGetLocalType()`, `PCASMType`, `PCCompositeType` 112612cd4985SMatthew G. Knepley @*/ 1127d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMSetLocalType(PC pc, PCCompositeType type) 1128d71ae5a4SJacob Faibussowitsch { 112912cd4985SMatthew G. Knepley PetscFunctionBegin; 113012cd4985SMatthew G. Knepley PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 113112cd4985SMatthew G. Knepley PetscValidLogicalCollectiveEnum(pc, type, 2); 1132cac4c232SBarry Smith PetscTryMethod(pc, "PCASMSetLocalType_C", (PC, PCCompositeType), (pc, type)); 11333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 113412cd4985SMatthew G. Knepley } 113512cd4985SMatthew G. Knepley 113612cd4985SMatthew G. Knepley /*@ 1137f1580f4eSBarry Smith PCASMGetLocalType - Gets the type of composition used for local problems in the additive Schwarz method, `PCASM`. 113812cd4985SMatthew G. Knepley 1139c3339decSBarry Smith Logically Collective 114012cd4985SMatthew G. Knepley 114112cd4985SMatthew G. Knepley Input Parameter: 114212cd4985SMatthew G. Knepley . pc - the preconditioner context 114312cd4985SMatthew G. Knepley 114412cd4985SMatthew G. Knepley Output Parameter: 114512cd4985SMatthew G. Knepley . type - type of composition, one of 114612cd4985SMatthew G. Knepley .vb 114712cd4985SMatthew G. Knepley PC_COMPOSITE_ADDITIVE - local additive combination 114812cd4985SMatthew G. Knepley PC_COMPOSITE_MULTIPLICATIVE - local multiplicative combination 114912cd4985SMatthew G. Knepley .ve 115012cd4985SMatthew G. Knepley 115112cd4985SMatthew G. Knepley Options Database Key: 115212cd4985SMatthew G. Knepley . -pc_asm_local_type [additive,multiplicative] - Sets local solver composition type 115312cd4985SMatthew G. Knepley 115412cd4985SMatthew G. Knepley Level: intermediate 115512cd4985SMatthew G. Knepley 1156562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMSetType()`, `PCASMGetType()`, `PCASMSetLocalType()`, `PCASMCreate()`, `PCASMType`, `PCCompositeType` 115712cd4985SMatthew G. Knepley @*/ 1158d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMGetLocalType(PC pc, PCCompositeType *type) 1159d71ae5a4SJacob Faibussowitsch { 116012cd4985SMatthew G. Knepley PetscFunctionBegin; 116112cd4985SMatthew G. Knepley PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 11624f572ea9SToby Isaac PetscAssertPointer(type, 2); 1163cac4c232SBarry Smith PetscUseMethod(pc, "PCASMGetLocalType_C", (PC, PCCompositeType *), (pc, type)); 11643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 116512cd4985SMatthew G. Knepley } 116612cd4985SMatthew G. Knepley 11676ed231c7SMatthew Knepley /*@ 11686ed231c7SMatthew Knepley PCASMSetSortIndices - Determines whether subdomain indices are sorted. 11696ed231c7SMatthew Knepley 1170c3339decSBarry Smith Logically Collective 11716ed231c7SMatthew Knepley 11726ed231c7SMatthew Knepley Input Parameters: 11736ed231c7SMatthew Knepley + pc - the preconditioner context 11746ed231c7SMatthew Knepley - doSort - sort the subdomain indices 11756ed231c7SMatthew Knepley 11766ed231c7SMatthew Knepley Level: intermediate 11776ed231c7SMatthew Knepley 1178562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMSetLocalSubdomains()`, `PCASMSetTotalSubdomains()`, `PCASMGetSubKSP()`, 1179db781477SPatrick Sanan `PCASMCreateSubdomains2D()` 11806ed231c7SMatthew Knepley @*/ 1181d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMSetSortIndices(PC pc, PetscBool doSort) 1182d71ae5a4SJacob Faibussowitsch { 11836ed231c7SMatthew Knepley PetscFunctionBegin; 11840700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1185acfcf0e5SJed Brown PetscValidLogicalCollectiveBool(pc, doSort, 2); 1186cac4c232SBarry Smith PetscTryMethod(pc, "PCASMSetSortIndices_C", (PC, PetscBool), (pc, doSort)); 11873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11886ed231c7SMatthew Knepley } 11896ed231c7SMatthew Knepley 11904b9ad928SBarry Smith /*@C 1191f1580f4eSBarry Smith PCASMGetSubKSP - Gets the local `KSP` contexts for all blocks on 11924b9ad928SBarry Smith this processor. 11934b9ad928SBarry Smith 1194c3339decSBarry Smith Collective iff first_local is requested 11954b9ad928SBarry Smith 11964b9ad928SBarry Smith Input Parameter: 11974b9ad928SBarry Smith . pc - the preconditioner context 11984b9ad928SBarry Smith 11994b9ad928SBarry Smith Output Parameters: 1200*a3b724e8SBarry Smith + n_local - the number of blocks on this processor or `NULL` 1201*a3b724e8SBarry Smith . first_local - the global number of the first block on this processor or `NULL`, all processors must request or all must pass `NULL` 1202f1580f4eSBarry Smith - ksp - the array of `KSP` contexts 12034b9ad928SBarry Smith 120420f4b53cSBarry Smith Level: advanced 120520f4b53cSBarry Smith 1206f1580f4eSBarry Smith Notes: 1207f1580f4eSBarry Smith After `PCASMGetSubKSP()` the array of `KSP`s is not to be freed. 12084b9ad928SBarry Smith 1209f1580f4eSBarry Smith You must call `KSPSetUp()` before calling `PCASMGetSubKSP()`. 12104b9ad928SBarry Smith 1211feefa0e1SJacob Faibussowitsch Fortran Notes: 1212f1580f4eSBarry Smith The output argument 'ksp' must be an array of sufficient length or `PETSC_NULL_KSP`. The latter can be used to learn the necessary length. 1213d29017ddSJed Brown 1214562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMSetTotalSubdomains()`, `PCASMSetOverlap()`, 1215db781477SPatrick Sanan `PCASMCreateSubdomains2D()`, 12164b9ad928SBarry Smith @*/ 1217d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMGetSubKSP(PC pc, PetscInt *n_local, PetscInt *first_local, KSP *ksp[]) 1218d71ae5a4SJacob Faibussowitsch { 12194b9ad928SBarry Smith PetscFunctionBegin; 12200700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1221cac4c232SBarry Smith PetscUseMethod(pc, "PCASMGetSubKSP_C", (PC, PetscInt *, PetscInt *, KSP **), (pc, n_local, first_local, ksp)); 12223ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 12234b9ad928SBarry Smith } 12244b9ad928SBarry Smith 12254b9ad928SBarry Smith /*MC 12263b09bd56SBarry Smith PCASM - Use the (restricted) additive Schwarz method, each block is (approximately) solved with 12271d27aa22SBarry Smith its own `KSP` object, {cite}`dryja1987additive` and {cite}`1sbg` 12284b9ad928SBarry Smith 12294b9ad928SBarry Smith Options Database Keys: 12301d27aa22SBarry Smith + -pc_asm_blocks <blks> - Sets total blocks. Defaults to one block per MPI process. 12314b9ad928SBarry Smith . -pc_asm_overlap <ovl> - Sets overlap 1232f1580f4eSBarry Smith . -pc_asm_type [basic,restrict,interpolate,none] - Sets `PCASMType`, default is restrict. See `PCASMSetType()` 123373ff1848SBarry Smith . -pc_asm_dm_subdomains <bool> - use subdomains defined by the `DM` with `DMCreateDomainDecomposition()` 1234f1580f4eSBarry Smith - -pc_asm_local_type [additive, multiplicative] - Sets `PCCompositeType`, default is additive. See `PCASMSetLocalType()` 12354b9ad928SBarry Smith 12364b9ad928SBarry Smith Level: beginner 12374b9ad928SBarry Smith 1238f1580f4eSBarry Smith Notes: 1239f1580f4eSBarry Smith If you run with, for example, 3 blocks on 1 processor or 3 blocks on 3 processors you 124073ff1848SBarry Smith will get a different convergence rate due to the default option of `-pc_asm_type restrict`. Use 124173ff1848SBarry Smith `-pc_asm_type basic` to get the same convergence behavior 1242f1580f4eSBarry Smith 1243f1580f4eSBarry Smith Each processor can have one or more blocks, but a block cannot be shared by more 1244f1580f4eSBarry Smith than one processor. Use `PCGASM` for subdomains shared by multiple processes. 1245f1580f4eSBarry Smith 124673ff1848SBarry Smith To set options on the solvers for each block append `-sub_` to all the `KSP`, and `PC` 124773ff1848SBarry Smith options database keys. For example, `-sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly` 1248f1580f4eSBarry Smith 1249f1580f4eSBarry Smith To set the options on the solvers separate for each block call `PCASMGetSubKSP()` 1250f1580f4eSBarry Smith and set the options directly on the resulting `KSP` object (you can access its `PC` with `KSPGetPC()`) 1251f1580f4eSBarry Smith 125273ff1848SBarry Smith If the `PC` has an associated `DM`, then, by default, `DMCreateDomainDecomposition()` is used to create the subdomains 125373ff1848SBarry Smith 1254562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCASMType`, `PCCompositeType`, 1255db781477SPatrick Sanan `PCBJACOBI`, `PCASMGetSubKSP()`, `PCASMSetLocalSubdomains()`, `PCASMType`, `PCASMGetType()`, `PCASMSetLocalType()`, `PCASMGetLocalType()` 1256db781477SPatrick Sanan `PCASMSetTotalSubdomains()`, `PCSetModifySubMatrices()`, `PCASMSetOverlap()`, `PCASMSetType()`, `PCCompositeType` 12574b9ad928SBarry Smith M*/ 12584b9ad928SBarry Smith 1259d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_ASM(PC pc) 1260d71ae5a4SJacob Faibussowitsch { 12614b9ad928SBarry Smith PC_ASM *osm; 12624b9ad928SBarry Smith 12634b9ad928SBarry Smith PetscFunctionBegin; 12644dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&osm)); 12652fa5cd67SKarl Rupp 12664b9ad928SBarry Smith osm->n = PETSC_DECIDE; 12674b9ad928SBarry Smith osm->n_local = 0; 12682b837212SDmitry Karpeev osm->n_local_true = PETSC_DECIDE; 12694b9ad928SBarry Smith osm->overlap = 1; 12700a545947SLisandro Dalcin osm->ksp = NULL; 12710a545947SLisandro Dalcin osm->restriction = NULL; 12720a545947SLisandro Dalcin osm->lprolongation = NULL; 12730a545947SLisandro Dalcin osm->lrestriction = NULL; 12740a545947SLisandro Dalcin osm->x = NULL; 12750a545947SLisandro Dalcin osm->y = NULL; 12760a545947SLisandro Dalcin osm->is = NULL; 12770a545947SLisandro Dalcin osm->is_local = NULL; 12780a545947SLisandro Dalcin osm->mat = NULL; 12790a545947SLisandro Dalcin osm->pmat = NULL; 12804b9ad928SBarry Smith osm->type = PC_ASM_RESTRICT; 128112cd4985SMatthew G. Knepley osm->loctype = PC_COMPOSITE_ADDITIVE; 12826ed231c7SMatthew Knepley osm->sort_indices = PETSC_TRUE; 1283d709ea83SDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 128480ec0b7dSPatrick Sanan osm->sub_mat_type = NULL; 12854b9ad928SBarry Smith 128692bb6962SLisandro Dalcin pc->data = (void *)osm; 12874b9ad928SBarry Smith pc->ops->apply = PCApply_ASM; 128848e38a8aSPierre Jolivet pc->ops->matapply = PCMatApply_ASM; 12894b9ad928SBarry Smith pc->ops->applytranspose = PCApplyTranspose_ASM; 12904b9ad928SBarry Smith pc->ops->setup = PCSetUp_ASM; 1291e91c6855SBarry Smith pc->ops->reset = PCReset_ASM; 12924b9ad928SBarry Smith pc->ops->destroy = PCDestroy_ASM; 12934b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_ASM; 12944b9ad928SBarry Smith pc->ops->setuponblocks = PCSetUpOnBlocks_ASM; 12954b9ad928SBarry Smith pc->ops->view = PCView_ASM; 12960a545947SLisandro Dalcin pc->ops->applyrichardson = NULL; 12974b9ad928SBarry Smith 12989566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetLocalSubdomains_C", PCASMSetLocalSubdomains_ASM)); 12999566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetTotalSubdomains_C", PCASMSetTotalSubdomains_ASM)); 13009566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetOverlap_C", PCASMSetOverlap_ASM)); 13019566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetType_C", PCASMSetType_ASM)); 13029566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMGetType_C", PCASMGetType_ASM)); 13039566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetLocalType_C", PCASMSetLocalType_ASM)); 13049566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMGetLocalType_C", PCASMGetLocalType_ASM)); 13059566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetSortIndices_C", PCASMSetSortIndices_ASM)); 13069566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMGetSubKSP_C", PCASMGetSubKSP_ASM)); 13079566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMGetSubMatType_C", PCASMGetSubMatType_ASM)); 13089566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetSubMatType_C", PCASMSetSubMatType_ASM)); 13093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13104b9ad928SBarry Smith } 13114b9ad928SBarry Smith 131292bb6962SLisandro Dalcin /*@C 131392bb6962SLisandro Dalcin PCASMCreateSubdomains - Creates the index sets for the overlapping Schwarz 1314f1580f4eSBarry Smith preconditioner, `PCASM`, for any problem on a general grid. 131592bb6962SLisandro Dalcin 131692bb6962SLisandro Dalcin Collective 131792bb6962SLisandro Dalcin 131892bb6962SLisandro Dalcin Input Parameters: 131992bb6962SLisandro Dalcin + A - The global matrix operator 132092bb6962SLisandro Dalcin - n - the number of local blocks 132192bb6962SLisandro Dalcin 13222fe279fdSBarry Smith Output Parameter: 132392bb6962SLisandro Dalcin . outis - the array of index sets defining the subdomains 132492bb6962SLisandro Dalcin 132592bb6962SLisandro Dalcin Level: advanced 132692bb6962SLisandro Dalcin 1327f1580f4eSBarry Smith Note: 1328f1580f4eSBarry Smith This generates nonoverlapping subdomains; the `PCASM` will generate the overlap 1329f1580f4eSBarry Smith from these if you use `PCASMSetLocalSubdomains()` 13307d6bfa3bSBarry Smith 1331feefa0e1SJacob Faibussowitsch Fortran Notes: 1332f1580f4eSBarry Smith You must provide the array outis[] already allocated of length n. 13337d6bfa3bSBarry Smith 1334562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMSetLocalSubdomains()`, `PCASMDestroySubdomains()` 133592bb6962SLisandro Dalcin @*/ 1336d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMCreateSubdomains(Mat A, PetscInt n, IS *outis[]) 1337d71ae5a4SJacob Faibussowitsch { 133892bb6962SLisandro Dalcin MatPartitioning mpart; 133992bb6962SLisandro Dalcin const char *prefix; 134092bb6962SLisandro Dalcin PetscInt i, j, rstart, rend, bs; 1341976e8c5aSLisandro Dalcin PetscBool hasop, isbaij = PETSC_FALSE, foundpart = PETSC_FALSE; 13420298fd71SBarry Smith Mat Ad = NULL, adj; 134392bb6962SLisandro Dalcin IS ispart, isnumb, *is; 134492bb6962SLisandro Dalcin 134592bb6962SLisandro Dalcin PetscFunctionBegin; 13460700a824SBarry Smith PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 13474f572ea9SToby Isaac PetscAssertPointer(outis, 3); 134863a3b9bcSJacob Faibussowitsch PetscCheck(n >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "number of local blocks must be > 0, n = %" PetscInt_FMT, n); 134992bb6962SLisandro Dalcin 135092bb6962SLisandro Dalcin /* Get prefix, row distribution, and block size */ 13519566063dSJacob Faibussowitsch PetscCall(MatGetOptionsPrefix(A, &prefix)); 13529566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A, &rstart, &rend)); 13539566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(A, &bs)); 135463a3b9bcSJacob 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); 135565e19b50SBarry Smith 135692bb6962SLisandro Dalcin /* Get diagonal block from matrix if possible */ 13579566063dSJacob Faibussowitsch PetscCall(MatHasOperation(A, MATOP_GET_DIAGONAL_BLOCK, &hasop)); 135848a46eb9SPierre Jolivet if (hasop) PetscCall(MatGetDiagonalBlock(A, &Ad)); 135992bb6962SLisandro Dalcin if (Ad) { 13609566063dSJacob Faibussowitsch PetscCall(PetscObjectBaseTypeCompare((PetscObject)Ad, MATSEQBAIJ, &isbaij)); 13619566063dSJacob Faibussowitsch if (!isbaij) PetscCall(PetscObjectBaseTypeCompare((PetscObject)Ad, MATSEQSBAIJ, &isbaij)); 136292bb6962SLisandro Dalcin } 136392bb6962SLisandro Dalcin if (Ad && n > 1) { 1364ace3abfcSBarry Smith PetscBool match, done; 136592bb6962SLisandro Dalcin /* Try to setup a good matrix partitioning if available */ 13669566063dSJacob Faibussowitsch PetscCall(MatPartitioningCreate(PETSC_COMM_SELF, &mpart)); 13679566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)mpart, prefix)); 13689566063dSJacob Faibussowitsch PetscCall(MatPartitioningSetFromOptions(mpart)); 13699566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)mpart, MATPARTITIONINGCURRENT, &match)); 137048a46eb9SPierre Jolivet if (!match) PetscCall(PetscObjectTypeCompare((PetscObject)mpart, MATPARTITIONINGSQUARE, &match)); 137192bb6962SLisandro Dalcin if (!match) { /* assume a "good" partitioner is available */ 13721a83f524SJed Brown PetscInt na; 13731a83f524SJed Brown const PetscInt *ia, *ja; 13749566063dSJacob Faibussowitsch PetscCall(MatGetRowIJ(Ad, 0, PETSC_TRUE, isbaij, &na, &ia, &ja, &done)); 137592bb6962SLisandro Dalcin if (done) { 137692bb6962SLisandro Dalcin /* Build adjacency matrix by hand. Unfortunately a call to 137792bb6962SLisandro Dalcin MatConvert(Ad,MATMPIADJ,MAT_INITIAL_MATRIX,&adj) will 137892bb6962SLisandro Dalcin remove the block-aij structure and we cannot expect 137992bb6962SLisandro Dalcin MatPartitioning to split vertices as we need */ 13800a545947SLisandro Dalcin PetscInt i, j, len, nnz, cnt, *iia = NULL, *jja = NULL; 13811a83f524SJed Brown const PetscInt *row; 138292bb6962SLisandro Dalcin nnz = 0; 138392bb6962SLisandro Dalcin for (i = 0; i < na; i++) { /* count number of nonzeros */ 138492bb6962SLisandro Dalcin len = ia[i + 1] - ia[i]; 138592bb6962SLisandro Dalcin row = ja + ia[i]; 138692bb6962SLisandro Dalcin for (j = 0; j < len; j++) { 138792bb6962SLisandro Dalcin if (row[j] == i) { /* don't count diagonal */ 13889371c9d4SSatish Balay len--; 13899371c9d4SSatish Balay break; 139092bb6962SLisandro Dalcin } 139192bb6962SLisandro Dalcin } 139292bb6962SLisandro Dalcin nnz += len; 139392bb6962SLisandro Dalcin } 13949566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(na + 1, &iia)); 13959566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nnz, &jja)); 139692bb6962SLisandro Dalcin nnz = 0; 139792bb6962SLisandro Dalcin iia[0] = 0; 139892bb6962SLisandro Dalcin for (i = 0; i < na; i++) { /* fill adjacency */ 139992bb6962SLisandro Dalcin cnt = 0; 140092bb6962SLisandro Dalcin len = ia[i + 1] - ia[i]; 140192bb6962SLisandro Dalcin row = ja + ia[i]; 140292bb6962SLisandro Dalcin for (j = 0; j < len; j++) { 140392bb6962SLisandro Dalcin if (row[j] != i) { /* if not diagonal */ 140492bb6962SLisandro Dalcin jja[nnz + cnt++] = row[j]; 140592bb6962SLisandro Dalcin } 140692bb6962SLisandro Dalcin } 140792bb6962SLisandro Dalcin nnz += cnt; 140892bb6962SLisandro Dalcin iia[i + 1] = nnz; 140992bb6962SLisandro Dalcin } 141092bb6962SLisandro Dalcin /* Partitioning of the adjacency matrix */ 14119566063dSJacob Faibussowitsch PetscCall(MatCreateMPIAdj(PETSC_COMM_SELF, na, na, iia, jja, NULL, &adj)); 14129566063dSJacob Faibussowitsch PetscCall(MatPartitioningSetAdjacency(mpart, adj)); 14139566063dSJacob Faibussowitsch PetscCall(MatPartitioningSetNParts(mpart, n)); 14149566063dSJacob Faibussowitsch PetscCall(MatPartitioningApply(mpart, &ispart)); 14159566063dSJacob Faibussowitsch PetscCall(ISPartitioningToNumbering(ispart, &isnumb)); 14169566063dSJacob Faibussowitsch PetscCall(MatDestroy(&adj)); 141792bb6962SLisandro Dalcin foundpart = PETSC_TRUE; 141892bb6962SLisandro Dalcin } 14199566063dSJacob Faibussowitsch PetscCall(MatRestoreRowIJ(Ad, 0, PETSC_TRUE, isbaij, &na, &ia, &ja, &done)); 142092bb6962SLisandro Dalcin } 14219566063dSJacob Faibussowitsch PetscCall(MatPartitioningDestroy(&mpart)); 142292bb6962SLisandro Dalcin } 142392bb6962SLisandro Dalcin 14249566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &is)); 142592bb6962SLisandro Dalcin *outis = is; 142692bb6962SLisandro Dalcin 142792bb6962SLisandro Dalcin if (!foundpart) { 142892bb6962SLisandro Dalcin /* Partitioning by contiguous chunks of rows */ 142992bb6962SLisandro Dalcin 143092bb6962SLisandro Dalcin PetscInt mbs = (rend - rstart) / bs; 143192bb6962SLisandro Dalcin PetscInt start = rstart; 143292bb6962SLisandro Dalcin for (i = 0; i < n; i++) { 143392bb6962SLisandro Dalcin PetscInt count = (mbs / n + ((mbs % n) > i)) * bs; 14349566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, count, start, 1, &is[i])); 143592bb6962SLisandro Dalcin start += count; 143692bb6962SLisandro Dalcin } 143792bb6962SLisandro Dalcin 143892bb6962SLisandro Dalcin } else { 143992bb6962SLisandro Dalcin /* Partitioning by adjacency of diagonal block */ 144092bb6962SLisandro Dalcin 144192bb6962SLisandro Dalcin const PetscInt *numbering; 144292bb6962SLisandro Dalcin PetscInt *count, nidx, *indices, *newidx, start = 0; 144392bb6962SLisandro Dalcin /* Get node count in each partition */ 14449566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &count)); 14459566063dSJacob Faibussowitsch PetscCall(ISPartitioningCount(ispart, n, count)); 144692bb6962SLisandro Dalcin if (isbaij && bs > 1) { /* adjust for the block-aij case */ 144792bb6962SLisandro Dalcin for (i = 0; i < n; i++) count[i] *= bs; 144892bb6962SLisandro Dalcin } 144992bb6962SLisandro Dalcin /* Build indices from node numbering */ 14509566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(isnumb, &nidx)); 14519566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nidx, &indices)); 145292bb6962SLisandro Dalcin for (i = 0; i < nidx; i++) indices[i] = i; /* needs to be initialized */ 14539566063dSJacob Faibussowitsch PetscCall(ISGetIndices(isnumb, &numbering)); 14549566063dSJacob Faibussowitsch PetscCall(PetscSortIntWithPermutation(nidx, numbering, indices)); 14559566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(isnumb, &numbering)); 145692bb6962SLisandro Dalcin if (isbaij && bs > 1) { /* adjust for the block-aij case */ 14579566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nidx * bs, &newidx)); 14582fa5cd67SKarl Rupp for (i = 0; i < nidx; i++) { 14592fa5cd67SKarl Rupp for (j = 0; j < bs; j++) newidx[i * bs + j] = indices[i] * bs + j; 14602fa5cd67SKarl Rupp } 14619566063dSJacob Faibussowitsch PetscCall(PetscFree(indices)); 146292bb6962SLisandro Dalcin nidx *= bs; 146392bb6962SLisandro Dalcin indices = newidx; 146492bb6962SLisandro Dalcin } 146592bb6962SLisandro Dalcin /* Shift to get global indices */ 146692bb6962SLisandro Dalcin for (i = 0; i < nidx; i++) indices[i] += rstart; 146792bb6962SLisandro Dalcin 146892bb6962SLisandro Dalcin /* Build the index sets for each block */ 146992bb6962SLisandro Dalcin for (i = 0; i < n; i++) { 14709566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, count[i], &indices[start], PETSC_COPY_VALUES, &is[i])); 14719566063dSJacob Faibussowitsch PetscCall(ISSort(is[i])); 147292bb6962SLisandro Dalcin start += count[i]; 147392bb6962SLisandro Dalcin } 147492bb6962SLisandro Dalcin 14759566063dSJacob Faibussowitsch PetscCall(PetscFree(count)); 14769566063dSJacob Faibussowitsch PetscCall(PetscFree(indices)); 14779566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isnumb)); 14789566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ispart)); 147992bb6962SLisandro Dalcin } 14803ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 148192bb6962SLisandro Dalcin } 148292bb6962SLisandro Dalcin 148392bb6962SLisandro Dalcin /*@C 148492bb6962SLisandro Dalcin PCASMDestroySubdomains - Destroys the index sets created with 1485f1580f4eSBarry Smith `PCASMCreateSubdomains()`. Should be called after setting subdomains with `PCASMSetLocalSubdomains()`. 148692bb6962SLisandro Dalcin 148792bb6962SLisandro Dalcin Collective 148892bb6962SLisandro Dalcin 148992bb6962SLisandro Dalcin Input Parameters: 149092bb6962SLisandro Dalcin + n - the number of index sets 14912b691e39SMatthew Knepley . is - the array of index sets 14922fe279fdSBarry Smith - is_local - the array of local index sets, can be `NULL` 149392bb6962SLisandro Dalcin 149492bb6962SLisandro Dalcin Level: advanced 149592bb6962SLisandro Dalcin 1496562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMCreateSubdomains()`, `PCASMSetLocalSubdomains()` 149792bb6962SLisandro Dalcin @*/ 1498d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMDestroySubdomains(PetscInt n, IS is[], IS is_local[]) 1499d71ae5a4SJacob Faibussowitsch { 150092bb6962SLisandro Dalcin PetscInt i; 15015fd66863SKarl Rupp 150292bb6962SLisandro Dalcin PetscFunctionBegin; 15033ba16761SJacob Faibussowitsch if (n <= 0) PetscFunctionReturn(PETSC_SUCCESS); 1504a35b7b57SDmitry Karpeev if (is) { 15054f572ea9SToby Isaac PetscAssertPointer(is, 2); 15069566063dSJacob Faibussowitsch for (i = 0; i < n; i++) PetscCall(ISDestroy(&is[i])); 15079566063dSJacob Faibussowitsch PetscCall(PetscFree(is)); 1508a35b7b57SDmitry Karpeev } 15092b691e39SMatthew Knepley if (is_local) { 15104f572ea9SToby Isaac PetscAssertPointer(is_local, 3); 15119566063dSJacob Faibussowitsch for (i = 0; i < n; i++) PetscCall(ISDestroy(&is_local[i])); 15129566063dSJacob Faibussowitsch PetscCall(PetscFree(is_local)); 15132b691e39SMatthew Knepley } 15143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 151592bb6962SLisandro Dalcin } 151692bb6962SLisandro Dalcin 15176141accfSBarry Smith /*@C 15184b9ad928SBarry Smith PCASMCreateSubdomains2D - Creates the index sets for the overlapping Schwarz 1519f1580f4eSBarry Smith preconditioner, `PCASM`, for a two-dimensional problem on a regular grid. 15204b9ad928SBarry Smith 15214b9ad928SBarry Smith Not Collective 15224b9ad928SBarry Smith 15234b9ad928SBarry Smith Input Parameters: 15246b867d5aSJose E. Roman + m - the number of mesh points in the x direction 15256b867d5aSJose E. Roman . n - the number of mesh points in the y direction 15266b867d5aSJose E. Roman . M - the number of subdomains in the x direction 15276b867d5aSJose E. Roman . N - the number of subdomains in the y direction 15284b9ad928SBarry Smith . dof - degrees of freedom per node 15294b9ad928SBarry Smith - overlap - overlap in mesh lines 15304b9ad928SBarry Smith 15314b9ad928SBarry Smith Output Parameters: 15324b9ad928SBarry Smith + Nsub - the number of subdomains created 15333d03bb48SJed Brown . is - array of index sets defining overlapping (if overlap > 0) subdomains 15343d03bb48SJed Brown - is_local - array of index sets defining non-overlapping subdomains 15354b9ad928SBarry Smith 153620f4b53cSBarry Smith Level: advanced 153720f4b53cSBarry Smith 15384b9ad928SBarry Smith Note: 1539f1580f4eSBarry Smith Presently `PCAMSCreateSubdomains2d()` is valid only for sequential 15404b9ad928SBarry Smith preconditioners. More general related routines are 1541f1580f4eSBarry Smith `PCASMSetTotalSubdomains()` and `PCASMSetLocalSubdomains()`. 15424b9ad928SBarry Smith 1543feefa0e1SJacob Faibussowitsch Fortran Notes: 15446141accfSBarry Smith The `IS` must be declared as an array of length long enough to hold `Nsub` entries 15456141accfSBarry Smith 1546562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMSetTotalSubdomains()`, `PCASMSetLocalSubdomains()`, `PCASMGetSubKSP()`, 1547db781477SPatrick Sanan `PCASMSetOverlap()` 15484b9ad928SBarry Smith @*/ 1549d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMCreateSubdomains2D(PetscInt m, PetscInt n, PetscInt M, PetscInt N, PetscInt dof, PetscInt overlap, PetscInt *Nsub, IS **is, IS **is_local) 1550d71ae5a4SJacob Faibussowitsch { 15513d03bb48SJed Brown PetscInt i, j, height, width, ystart, xstart, yleft, yright, xleft, xright, loc_outer; 155213f74950SBarry Smith PetscInt nidx, *idx, loc, ii, jj, count; 15534b9ad928SBarry Smith 15544b9ad928SBarry Smith PetscFunctionBegin; 15557827d75bSBarry Smith PetscCheck(dof == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "dof must be 1"); 15564b9ad928SBarry Smith 15574b9ad928SBarry Smith *Nsub = N * M; 15589566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(*Nsub, is)); 15599566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(*Nsub, is_local)); 15604b9ad928SBarry Smith ystart = 0; 15613d03bb48SJed Brown loc_outer = 0; 15624b9ad928SBarry Smith for (i = 0; i < N; i++) { 15634b9ad928SBarry Smith height = n / N + ((n % N) > i); /* height of subdomain */ 15647827d75bSBarry Smith PetscCheck(height >= 2, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many N subdomains for mesh dimension n"); 15659371c9d4SSatish Balay yleft = ystart - overlap; 15669371c9d4SSatish Balay if (yleft < 0) yleft = 0; 15679371c9d4SSatish Balay yright = ystart + height + overlap; 15689371c9d4SSatish Balay if (yright > n) yright = n; 15694b9ad928SBarry Smith xstart = 0; 15704b9ad928SBarry Smith for (j = 0; j < M; j++) { 15714b9ad928SBarry Smith width = m / M + ((m % M) > j); /* width of subdomain */ 15727827d75bSBarry Smith PetscCheck(width >= 2, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many M subdomains for mesh dimension m"); 15739371c9d4SSatish Balay xleft = xstart - overlap; 15749371c9d4SSatish Balay if (xleft < 0) xleft = 0; 15759371c9d4SSatish Balay xright = xstart + width + overlap; 15769371c9d4SSatish Balay if (xright > m) xright = m; 15774b9ad928SBarry Smith nidx = (xright - xleft) * (yright - yleft); 15789566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nidx, &idx)); 15794b9ad928SBarry Smith loc = 0; 15804b9ad928SBarry Smith for (ii = yleft; ii < yright; ii++) { 15814b9ad928SBarry Smith count = m * ii + xleft; 15822fa5cd67SKarl Rupp for (jj = xleft; jj < xright; jj++) idx[loc++] = count++; 15834b9ad928SBarry Smith } 15849566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, nidx, idx, PETSC_COPY_VALUES, (*is) + loc_outer)); 15853d03bb48SJed Brown if (overlap == 0) { 15869566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)(*is)[loc_outer])); 15872fa5cd67SKarl Rupp 15883d03bb48SJed Brown (*is_local)[loc_outer] = (*is)[loc_outer]; 15893d03bb48SJed Brown } else { 15903d03bb48SJed Brown for (loc = 0, ii = ystart; ii < ystart + height; ii++) { 1591ad540459SPierre Jolivet for (jj = xstart; jj < xstart + width; jj++) idx[loc++] = m * ii + jj; 15923d03bb48SJed Brown } 15939566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, loc, idx, PETSC_COPY_VALUES, *is_local + loc_outer)); 15943d03bb48SJed Brown } 15959566063dSJacob Faibussowitsch PetscCall(PetscFree(idx)); 15964b9ad928SBarry Smith xstart += width; 15973d03bb48SJed Brown loc_outer++; 15984b9ad928SBarry Smith } 15994b9ad928SBarry Smith ystart += height; 16004b9ad928SBarry Smith } 16019566063dSJacob Faibussowitsch for (i = 0; i < *Nsub; i++) PetscCall(ISSort((*is)[i])); 16023ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 16034b9ad928SBarry Smith } 16044b9ad928SBarry Smith 16054b9ad928SBarry Smith /*@C 16064b9ad928SBarry Smith PCASMGetLocalSubdomains - Gets the local subdomains (for this processor 1607f1580f4eSBarry Smith only) for the additive Schwarz preconditioner, `PCASM`. 16084b9ad928SBarry Smith 1609ad4df100SBarry Smith Not Collective 16104b9ad928SBarry Smith 16114b9ad928SBarry Smith Input Parameter: 16124b9ad928SBarry Smith . pc - the preconditioner context 16134b9ad928SBarry Smith 16144b9ad928SBarry Smith Output Parameters: 1615f17a6784SPierre Jolivet + n - if requested, the number of subdomains for this processor (default value = 1) 1616f17a6784SPierre Jolivet . is - if requested, the index sets that define the subdomains for this processor 161720f4b53cSBarry Smith - is_local - if requested, the index sets that define the local part of the subdomains for this processor (can be `NULL`) 161820f4b53cSBarry Smith 161920f4b53cSBarry Smith Level: advanced 16204b9ad928SBarry Smith 1621f1580f4eSBarry Smith Note: 1622f1580f4eSBarry Smith The `IS` numbering is in the parallel, global numbering of the vector. 16234b9ad928SBarry Smith 1624562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMSetTotalSubdomains()`, `PCASMSetOverlap()`, `PCASMGetSubKSP()`, 1625db781477SPatrick Sanan `PCASMCreateSubdomains2D()`, `PCASMSetLocalSubdomains()`, `PCASMGetLocalSubmatrices()` 16264b9ad928SBarry Smith @*/ 1627d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMGetLocalSubdomains(PC pc, PetscInt *n, IS *is[], IS *is_local[]) 1628d71ae5a4SJacob Faibussowitsch { 16292a808120SBarry Smith PC_ASM *osm = (PC_ASM *)pc->data; 1630ace3abfcSBarry Smith PetscBool match; 16314b9ad928SBarry Smith 16324b9ad928SBarry Smith PetscFunctionBegin; 16330700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 16344f572ea9SToby Isaac if (n) PetscAssertPointer(n, 2); 16354f572ea9SToby Isaac if (is) PetscAssertPointer(is, 3); 16364f572ea9SToby Isaac if (is_local) PetscAssertPointer(is_local, 4); 16379566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCASM, &match)); 163828b400f6SJacob Faibussowitsch PetscCheck(match, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "PC is not a PCASM"); 16394b9ad928SBarry Smith if (n) *n = osm->n_local_true; 16404b9ad928SBarry Smith if (is) *is = osm->is; 16412b691e39SMatthew Knepley if (is_local) *is_local = osm->is_local; 16423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 16434b9ad928SBarry Smith } 16444b9ad928SBarry Smith 16454b9ad928SBarry Smith /*@C 16464b9ad928SBarry Smith PCASMGetLocalSubmatrices - Gets the local submatrices (for this processor 1647f1580f4eSBarry Smith only) for the additive Schwarz preconditioner, `PCASM`. 16484b9ad928SBarry Smith 1649ad4df100SBarry Smith Not Collective 16504b9ad928SBarry Smith 16514b9ad928SBarry Smith Input Parameter: 16524b9ad928SBarry Smith . pc - the preconditioner context 16534b9ad928SBarry Smith 16544b9ad928SBarry Smith Output Parameters: 1655f17a6784SPierre Jolivet + n - if requested, the number of matrices for this processor (default value = 1) 1656f17a6784SPierre Jolivet - mat - if requested, the matrices 16574b9ad928SBarry Smith 16584b9ad928SBarry Smith Level: advanced 16594b9ad928SBarry Smith 166095452b02SPatrick Sanan Notes: 1661f1580f4eSBarry Smith Call after `PCSetUp()` (or `KSPSetUp()`) but before `PCApply()` and before `PCSetUpOnBlocks()`) 1662cf739d55SBarry Smith 1663f1580f4eSBarry Smith Usually one would use `PCSetModifySubMatrices()` to change the submatrices in building the preconditioner. 1664cf739d55SBarry Smith 1665562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMSetTotalSubdomains()`, `PCASMSetOverlap()`, `PCASMGetSubKSP()`, 1666db781477SPatrick Sanan `PCASMCreateSubdomains2D()`, `PCASMSetLocalSubdomains()`, `PCASMGetLocalSubdomains()`, `PCSetModifySubMatrices()` 16674b9ad928SBarry Smith @*/ 1668d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMGetLocalSubmatrices(PC pc, PetscInt *n, Mat *mat[]) 1669d71ae5a4SJacob Faibussowitsch { 16704b9ad928SBarry Smith PC_ASM *osm; 1671ace3abfcSBarry Smith PetscBool match; 16724b9ad928SBarry Smith 16734b9ad928SBarry Smith PetscFunctionBegin; 16740700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 16754f572ea9SToby Isaac if (n) PetscAssertPointer(n, 2); 16764f572ea9SToby Isaac if (mat) PetscAssertPointer(mat, 3); 167728b400f6SJacob Faibussowitsch PetscCheck(pc->setupcalled, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must call after KSPSetUp() or PCSetUp()."); 16789566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCASM, &match)); 167992bb6962SLisandro Dalcin if (!match) { 168092bb6962SLisandro Dalcin if (n) *n = 0; 16810298fd71SBarry Smith if (mat) *mat = NULL; 168292bb6962SLisandro Dalcin } else { 16834b9ad928SBarry Smith osm = (PC_ASM *)pc->data; 16844b9ad928SBarry Smith if (n) *n = osm->n_local_true; 16854b9ad928SBarry Smith if (mat) *mat = osm->pmat; 168692bb6962SLisandro Dalcin } 16873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 16884b9ad928SBarry Smith } 1689d709ea83SDmitry Karpeev 1690d709ea83SDmitry Karpeev /*@ 1691f1580f4eSBarry Smith PCASMSetDMSubdomains - Indicates whether to use `DMCreateDomainDecomposition()` to define the subdomains, whenever possible. 1692f1ee410cSBarry Smith 1693d709ea83SDmitry Karpeev Logically Collective 1694d709ea83SDmitry Karpeev 1695d8d19677SJose E. Roman Input Parameters: 1696d709ea83SDmitry Karpeev + pc - the preconditioner 1697f1580f4eSBarry Smith - flg - boolean indicating whether to use subdomains defined by the `DM` 1698d709ea83SDmitry Karpeev 1699d709ea83SDmitry Karpeev Options Database Key: 170073ff1848SBarry Smith . -pc_asm_dm_subdomains <bool> - use subdomains defined by the `DM` with `DMCreateDomainDecomposition()` 1701d709ea83SDmitry Karpeev 1702d709ea83SDmitry Karpeev Level: intermediate 1703d709ea83SDmitry Karpeev 1704f1580f4eSBarry Smith Note: 1705f1580f4eSBarry Smith `PCASMSetTotalSubdomains()` and `PCASMSetOverlap()` take precedence over `PCASMSetDMSubdomains()`, 1706d709ea83SDmitry Karpeev so setting either of the first two effectively turns the latter off. 1707d709ea83SDmitry Karpeev 170873ff1848SBarry Smith Developer Note: 170973ff1848SBarry Smith This should be `PCASMSetUseDMSubdomains()`, similarly for the options database key 171073ff1848SBarry Smith 1711562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMGetDMSubdomains()`, `PCASMSetTotalSubdomains()`, `PCASMSetOverlap()` 1712db781477SPatrick Sanan `PCASMCreateSubdomains2D()`, `PCASMSetLocalSubdomains()`, `PCASMGetLocalSubdomains()` 1713d709ea83SDmitry Karpeev @*/ 1714d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMSetDMSubdomains(PC pc, PetscBool flg) 1715d71ae5a4SJacob Faibussowitsch { 1716d709ea83SDmitry Karpeev PC_ASM *osm = (PC_ASM *)pc->data; 1717d709ea83SDmitry Karpeev PetscBool match; 1718d709ea83SDmitry Karpeev 1719d709ea83SDmitry Karpeev PetscFunctionBegin; 1720d709ea83SDmitry Karpeev PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1721d709ea83SDmitry Karpeev PetscValidLogicalCollectiveBool(pc, flg, 2); 172228b400f6SJacob Faibussowitsch PetscCheck(!pc->setupcalled, ((PetscObject)pc)->comm, PETSC_ERR_ARG_WRONGSTATE, "Not for a setup PC."); 17239566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCASM, &match)); 1724ad540459SPierre Jolivet if (match) osm->dm_subdomains = flg; 17253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1726d709ea83SDmitry Karpeev } 1727d709ea83SDmitry Karpeev 1728d709ea83SDmitry Karpeev /*@ 1729f1580f4eSBarry Smith PCASMGetDMSubdomains - Returns flag indicating whether to use `DMCreateDomainDecomposition()` to define the subdomains, whenever possible. 1730f1580f4eSBarry Smith 1731d709ea83SDmitry Karpeev Not Collective 1732d709ea83SDmitry Karpeev 1733d709ea83SDmitry Karpeev Input Parameter: 1734d709ea83SDmitry Karpeev . pc - the preconditioner 1735d709ea83SDmitry Karpeev 1736d709ea83SDmitry Karpeev Output Parameter: 173720f4b53cSBarry Smith . flg - boolean indicating whether to use subdomains defined by the `DM` 1738d709ea83SDmitry Karpeev 1739d709ea83SDmitry Karpeev Level: intermediate 1740d709ea83SDmitry Karpeev 174173ff1848SBarry Smith Developer Note: 174273ff1848SBarry Smith This should be `PCASMSetUseDMSubdomains()` 174373ff1848SBarry Smith 1744562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMSetDMSubdomains()`, `PCASMSetTotalSubdomains()`, `PCASMSetOverlap()` 1745db781477SPatrick Sanan `PCASMCreateSubdomains2D()`, `PCASMSetLocalSubdomains()`, `PCASMGetLocalSubdomains()` 1746d709ea83SDmitry Karpeev @*/ 1747d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMGetDMSubdomains(PC pc, PetscBool *flg) 1748d71ae5a4SJacob Faibussowitsch { 1749d709ea83SDmitry Karpeev PC_ASM *osm = (PC_ASM *)pc->data; 1750d709ea83SDmitry Karpeev PetscBool match; 1751d709ea83SDmitry Karpeev 1752d709ea83SDmitry Karpeev PetscFunctionBegin; 1753d709ea83SDmitry Karpeev PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 17544f572ea9SToby Isaac PetscAssertPointer(flg, 2); 17559566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCASM, &match)); 175656ea09ceSStefano Zampini if (match) *flg = osm->dm_subdomains; 175756ea09ceSStefano Zampini else *flg = PETSC_FALSE; 17583ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1759d709ea83SDmitry Karpeev } 176080ec0b7dSPatrick Sanan 1761feefa0e1SJacob Faibussowitsch /*@C 1762f1580f4eSBarry Smith PCASMGetSubMatType - Gets the matrix type used for `PCASM` subsolves, as a string. 176380ec0b7dSPatrick Sanan 176480ec0b7dSPatrick Sanan Not Collective 176580ec0b7dSPatrick Sanan 176680ec0b7dSPatrick Sanan Input Parameter: 1767f1580f4eSBarry Smith . pc - the `PC` 176880ec0b7dSPatrick Sanan 176980ec0b7dSPatrick Sanan Output Parameter: 1770feefa0e1SJacob Faibussowitsch . sub_mat_type - name of matrix type 177180ec0b7dSPatrick Sanan 177280ec0b7dSPatrick Sanan Level: advanced 177380ec0b7dSPatrick Sanan 1774562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMSetSubMatType()`, `PCSetType()`, `VecSetType()`, `MatType`, `Mat` 177580ec0b7dSPatrick Sanan @*/ 1776d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMGetSubMatType(PC pc, MatType *sub_mat_type) 1777d71ae5a4SJacob Faibussowitsch { 177856ea09ceSStefano Zampini PetscFunctionBegin; 177956ea09ceSStefano Zampini PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1780cac4c232SBarry Smith PetscTryMethod(pc, "PCASMGetSubMatType_C", (PC, MatType *), (pc, sub_mat_type)); 17813ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 178280ec0b7dSPatrick Sanan } 178380ec0b7dSPatrick Sanan 1784feefa0e1SJacob Faibussowitsch /*@C 1785f1580f4eSBarry Smith PCASMSetSubMatType - Set the type of matrix used for `PCASM` subsolves 178680ec0b7dSPatrick Sanan 1787c3339decSBarry Smith Collective 178880ec0b7dSPatrick Sanan 178980ec0b7dSPatrick Sanan Input Parameters: 1790f1580f4eSBarry Smith + pc - the `PC` object 1791f1580f4eSBarry Smith - sub_mat_type - the `MatType` 179280ec0b7dSPatrick Sanan 179380ec0b7dSPatrick Sanan Options Database Key: 1794f1580f4eSBarry Smith . -pc_asm_sub_mat_type <sub_mat_type> - Sets the matrix type used for subsolves, for example, seqaijviennacl. 1795f1580f4eSBarry Smith If you specify a base name like aijviennacl, the corresponding sequential type is assumed. 179680ec0b7dSPatrick Sanan 1797f1580f4eSBarry Smith Note: 17982fe279fdSBarry Smith See `MatType` for available types 179980ec0b7dSPatrick Sanan 180080ec0b7dSPatrick Sanan Level: advanced 180180ec0b7dSPatrick Sanan 1802562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMGetSubMatType()`, `PCSetType()`, `VecSetType()`, `MatType`, `Mat` 180380ec0b7dSPatrick Sanan @*/ 1804d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMSetSubMatType(PC pc, MatType sub_mat_type) 1805d71ae5a4SJacob Faibussowitsch { 180656ea09ceSStefano Zampini PetscFunctionBegin; 180756ea09ceSStefano Zampini PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1808cac4c232SBarry Smith PetscTryMethod(pc, "PCASMSetSubMatType_C", (PC, MatType), (pc, sub_mat_type)); 18093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 181080ec0b7dSPatrick Sanan } 1811