14b9ad928SBarry Smith /* 24b9ad928SBarry Smith This file defines an additive Schwarz preconditioner for any Mat implementation. 34b9ad928SBarry Smith 44b9ad928SBarry Smith Note that each processor may have any number of subdomains. But in order to 54b9ad928SBarry Smith deal easily with the VecScatter(), we treat each processor as if it has the 64b9ad928SBarry Smith same number of subdomains. 74b9ad928SBarry Smith 84b9ad928SBarry Smith n - total number of true subdomains on all processors 94b9ad928SBarry Smith n_local_true - actual number of subdomains on this processor 104b9ad928SBarry Smith n_local = maximum over all processors of n_local_true 114b9ad928SBarry Smith */ 124b9ad928SBarry Smith 1334eca32bSPierre Jolivet #include <petsc/private/pcasmimpl.h> /*I "petscpc.h" I*/ 1426cc229bSBarry Smith #include "petsc/private/matimpl.h" 154b9ad928SBarry Smith 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) { 469566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 479566063dSJacob Faibussowitsch PetscCall(KSPView(osm->ksp[0], sviewer)); 489566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 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)); 6263a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(sviewer, "[%d] local block number %" PetscInt_FMT ", size = %" PetscInt_FMT "\n", (int)rank, i, bsz)); 639566063dSJacob Faibussowitsch PetscCall(KSPView(osm->ksp[i], sviewer)); 649566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(sviewer, "- - - - - - - - - - - - - - - - - -\n")); 654b9ad928SBarry Smith } 669566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer)); 679566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 689566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 699566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 704b9ad928SBarry Smith } 714b9ad928SBarry Smith } else if (isstring) { 7263a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerStringSPrintf(viewer, " blocks=%" PetscInt_FMT ", overlap=%" PetscInt_FMT ", type=%s", osm->n, osm->overlap, PCASMTypes[osm->type])); 739566063dSJacob Faibussowitsch PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer)); 749566063dSJacob Faibussowitsch if (osm->ksp) PetscCall(KSPView(osm->ksp[0], sviewer)); 759566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer)); 764b9ad928SBarry Smith } 773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 784b9ad928SBarry Smith } 794b9ad928SBarry Smith 80d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCASMPrintSubdomains(PC pc) 81d71ae5a4SJacob Faibussowitsch { 8292bb6962SLisandro Dalcin PC_ASM *osm = (PC_ASM *)pc->data; 8392bb6962SLisandro Dalcin const char *prefix; 8492bb6962SLisandro Dalcin char fname[PETSC_MAX_PATH_LEN + 1]; 85643535aeSDmitry Karpeev PetscViewer viewer, sviewer; 86643535aeSDmitry Karpeev char *s; 8792bb6962SLisandro Dalcin PetscInt i, j, nidx; 8892bb6962SLisandro Dalcin const PetscInt *idx; 89643535aeSDmitry Karpeev PetscMPIInt rank, size; 90846783a0SBarry Smith 9192bb6962SLisandro Dalcin PetscFunctionBegin; 929566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size)); 939566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank)); 949566063dSJacob Faibussowitsch PetscCall(PCGetOptionsPrefix(pc, &prefix)); 959566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, prefix, "-pc_asm_print_subdomains", fname, sizeof(fname), NULL)); 96c6a7a370SJeremy L Thompson if (fname[0] == 0) PetscCall(PetscStrncpy(fname, "stdout", sizeof(fname))); 979566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIOpen(PetscObjectComm((PetscObject)pc), fname, &viewer)); 98643535aeSDmitry Karpeev for (i = 0; i < osm->n_local; i++) { 99643535aeSDmitry Karpeev if (i < osm->n_local_true) { 1009566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(osm->is[i], &nidx)); 1019566063dSJacob Faibussowitsch PetscCall(ISGetIndices(osm->is[i], &idx)); 102643535aeSDmitry Karpeev /* Print to a string viewer; no more than 15 characters per index plus 512 char for the header.*/ 10336a9e3b9SBarry Smith #define len 16 * (nidx + 1) + 512 1049566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(len, &s)); 1059566063dSJacob Faibussowitsch PetscCall(PetscViewerStringOpen(PETSC_COMM_SELF, s, len, &sviewer)); 10636a9e3b9SBarry Smith #undef len 10763a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerStringSPrintf(sviewer, "[%d:%d] Subdomain %" PetscInt_FMT " with overlap:\n", rank, size, i)); 10848a46eb9SPierre Jolivet for (j = 0; j < nidx; j++) PetscCall(PetscViewerStringSPrintf(sviewer, "%" PetscInt_FMT " ", idx[j])); 1099566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(osm->is[i], &idx)); 1109566063dSJacob Faibussowitsch PetscCall(PetscViewerStringSPrintf(sviewer, "\n")); 1119566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&sviewer)); 1129566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 11363a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "%s", s)); 1149566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 1159566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 1169566063dSJacob Faibussowitsch PetscCall(PetscFree(s)); 1172b691e39SMatthew Knepley if (osm->is_local) { 118643535aeSDmitry Karpeev /* Print to a string viewer; no more than 15 characters per index plus 512 char for the header.*/ 11936a9e3b9SBarry Smith #define len 16 * (nidx + 1) + 512 1209566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(len, &s)); 1219566063dSJacob Faibussowitsch PetscCall(PetscViewerStringOpen(PETSC_COMM_SELF, s, len, &sviewer)); 12236a9e3b9SBarry Smith #undef len 12363a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerStringSPrintf(sviewer, "[%d:%d] Subdomain %" PetscInt_FMT " without overlap:\n", rank, size, i)); 1249566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(osm->is_local[i], &nidx)); 1259566063dSJacob Faibussowitsch PetscCall(ISGetIndices(osm->is_local[i], &idx)); 12648a46eb9SPierre Jolivet for (j = 0; j < nidx; j++) PetscCall(PetscViewerStringSPrintf(sviewer, "%" PetscInt_FMT " ", idx[j])); 1279566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(osm->is_local[i], &idx)); 1289566063dSJacob Faibussowitsch PetscCall(PetscViewerStringSPrintf(sviewer, "\n")); 1299566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&sviewer)); 1309566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 13163a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "%s", s)); 1329566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 1339566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 1349566063dSJacob Faibussowitsch PetscCall(PetscFree(s)); 135643535aeSDmitry Karpeev } 1362fa5cd67SKarl Rupp } else { 137643535aeSDmitry Karpeev /* Participate in collective viewer calls. */ 1389566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 1399566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 1409566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 141643535aeSDmitry Karpeev /* Assume either all ranks have is_local or none do. */ 142643535aeSDmitry Karpeev if (osm->is_local) { 1439566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 1449566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 1459566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 146643535aeSDmitry Karpeev } 147fdc48646SDmitry Karpeev } 14892bb6962SLisandro Dalcin } 1499566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 1509566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 1513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 15292bb6962SLisandro Dalcin } 15392bb6962SLisandro Dalcin 154d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUp_ASM(PC pc) 155d71ae5a4SJacob Faibussowitsch { 1564b9ad928SBarry Smith PC_ASM *osm = (PC_ASM *)pc->data; 15757501b6eSBarry Smith PetscBool flg; 15887e86712SBarry Smith PetscInt i, m, m_local; 1594b9ad928SBarry Smith MatReuse scall = MAT_REUSE_MATRIX; 1604b9ad928SBarry Smith IS isl; 1614b9ad928SBarry Smith KSP ksp; 1624b9ad928SBarry Smith PC subpc; 1632dcb1b2aSMatthew Knepley const char *prefix, *pprefix; 16423ce1328SBarry Smith Vec vec; 1650298fd71SBarry Smith DM *domain_dm = 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) { 2809566063dSJacob Faibussowitsch PetscCall(MatDestroyMatrices(osm->n_local_true, &osm->pmat)); 2814b9ad928SBarry Smith scall = MAT_INITIAL_MATRIX; 2824b9ad928SBarry Smith } 2834b9ad928SBarry Smith } 2844b9ad928SBarry Smith 285b58cb649SBarry Smith /* Destroy previous submatrices of a different type than pc->pmat since MAT_REUSE_MATRIX won't work in that case */ 28682b5ce2aSStefano Zampini if (scall == MAT_REUSE_MATRIX && osm->sub_mat_type) { 28748a46eb9SPierre Jolivet if (osm->n_local_true > 0) PetscCall(MatDestroySubMatrices(osm->n_local_true, &osm->pmat)); 288b58cb649SBarry Smith scall = MAT_INITIAL_MATRIX; 289b58cb649SBarry Smith } 290b58cb649SBarry Smith 29192bb6962SLisandro Dalcin /* 29292bb6962SLisandro Dalcin Extract out the submatrices 29392bb6962SLisandro Dalcin */ 2949566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(pc->pmat, osm->n_local_true, osm->is, osm->is, scall, &osm->pmat)); 29592bb6962SLisandro Dalcin if (scall == MAT_INITIAL_MATRIX) { 2969566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptionsPrefix((PetscObject)pc->pmat, &pprefix)); 297aa624791SPierre Jolivet for (i = 0; i < osm->n_local_true; i++) PetscCall(PetscObjectSetOptionsPrefix((PetscObject)osm->pmat[i], pprefix)); 29892bb6962SLisandro Dalcin } 29980ec0b7dSPatrick Sanan 30080ec0b7dSPatrick Sanan /* Convert the types of the submatrices (if needbe) */ 30180ec0b7dSPatrick Sanan if (osm->sub_mat_type) { 30248a46eb9SPierre Jolivet for (i = 0; i < osm->n_local_true; i++) PetscCall(MatConvert(osm->pmat[i], osm->sub_mat_type, MAT_INPLACE_MATRIX, &(osm->pmat[i]))); 30380ec0b7dSPatrick Sanan } 30480ec0b7dSPatrick Sanan 30580ec0b7dSPatrick Sanan if (!pc->setupcalled) { 30656ea09ceSStefano Zampini VecType vtype; 30756ea09ceSStefano Zampini 30880ec0b7dSPatrick Sanan /* Create the local work vectors (from the local matrices) and scatter contexts */ 3099566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(pc->pmat, &vec, NULL)); 3105b3c0d42Seaulisa 311677c7726SPierre 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"); 312677c7726SPierre Jolivet if (osm->is_local && osm->type != PC_ASM_BASIC && osm->loctype == PC_COMPOSITE_ADDITIVE) PetscCall(PetscMalloc1(osm->n_local_true, &osm->lprolongation)); 3139566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(osm->n_local_true, &osm->lrestriction)); 3149566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(osm->n_local_true, &osm->x)); 3159566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(osm->n_local_true, &osm->y)); 316347574c9Seaulisa 3179566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(osm->lis, &m)); 3189566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, m, 0, 1, &isl)); 3199566063dSJacob Faibussowitsch PetscCall(MatGetVecType(osm->pmat[0], &vtype)); 3209566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_SELF, &osm->lx)); 3219566063dSJacob Faibussowitsch PetscCall(VecSetSizes(osm->lx, m, m)); 3229566063dSJacob Faibussowitsch PetscCall(VecSetType(osm->lx, vtype)); 3239566063dSJacob Faibussowitsch PetscCall(VecDuplicate(osm->lx, &osm->ly)); 3249566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(vec, osm->lis, osm->lx, isl, &osm->restriction)); 3259566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isl)); 326347574c9Seaulisa 32780ec0b7dSPatrick Sanan for (i = 0; i < osm->n_local_true; ++i) { 3285b3c0d42Seaulisa ISLocalToGlobalMapping ltog; 3295b3c0d42Seaulisa IS isll; 3305b3c0d42Seaulisa const PetscInt *idx_is; 3315b3c0d42Seaulisa PetscInt *idx_lis, nout; 3325b3c0d42Seaulisa 3339566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(osm->is[i], &m)); 3349566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(osm->pmat[i], &osm->x[i], NULL)); 3359566063dSJacob Faibussowitsch PetscCall(VecDuplicate(osm->x[i], &osm->y[i])); 33655817e79SBarry Smith 337b0de9f23SBarry Smith /* generate a scatter from ly to y[i] picking all the overlapping is[i] entries */ 3389566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreateIS(osm->lis, <og)); 3399566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(osm->is[i], &m)); 3409566063dSJacob Faibussowitsch PetscCall(ISGetIndices(osm->is[i], &idx_is)); 3419566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m, &idx_lis)); 3429566063dSJacob Faibussowitsch PetscCall(ISGlobalToLocalMappingApply(ltog, IS_GTOLM_DROP, m, idx_is, &nout, idx_lis)); 3437827d75bSBarry Smith PetscCheck(nout == m, PETSC_COMM_SELF, PETSC_ERR_PLIB, "is not a subset of lis"); 3449566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(osm->is[i], &idx_is)); 3459566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, m, idx_lis, PETSC_OWN_POINTER, &isll)); 3469566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingDestroy(<og)); 3479566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, m, 0, 1, &isl)); 3489566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(osm->ly, isll, osm->y[i], isl, &osm->lrestriction[i])); 3499566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isll)); 3509566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isl)); 351910cf402Sprj- if (osm->lprolongation) { /* generate a scatter from y[i] to ly picking only the the non-overlapping is_local[i] entries */ 352d8b3b5e3Seaulisa ISLocalToGlobalMapping ltog; 353d8b3b5e3Seaulisa IS isll, isll_local; 354d8b3b5e3Seaulisa const PetscInt *idx_local; 355d8b3b5e3Seaulisa PetscInt *idx1, *idx2, nout; 356d8b3b5e3Seaulisa 3579566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(osm->is_local[i], &m_local)); 3589566063dSJacob Faibussowitsch PetscCall(ISGetIndices(osm->is_local[i], &idx_local)); 359d8b3b5e3Seaulisa 3609566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreateIS(osm->is[i], <og)); 3619566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m_local, &idx1)); 3629566063dSJacob Faibussowitsch PetscCall(ISGlobalToLocalMappingApply(ltog, IS_GTOLM_DROP, m_local, idx_local, &nout, idx1)); 3639566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingDestroy(<og)); 3647827d75bSBarry Smith PetscCheck(nout == m_local, PETSC_COMM_SELF, PETSC_ERR_PLIB, "is_local not a subset of is"); 3659566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, m_local, idx1, PETSC_OWN_POINTER, &isll)); 366d8b3b5e3Seaulisa 3679566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreateIS(osm->lis, <og)); 3689566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m_local, &idx2)); 3699566063dSJacob Faibussowitsch PetscCall(ISGlobalToLocalMappingApply(ltog, IS_GTOLM_DROP, m_local, idx_local, &nout, idx2)); 3709566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingDestroy(<og)); 3717827d75bSBarry Smith PetscCheck(nout == m_local, PETSC_COMM_SELF, PETSC_ERR_PLIB, "is_local not a subset of lis"); 3729566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, m_local, idx2, PETSC_OWN_POINTER, &isll_local)); 373d8b3b5e3Seaulisa 3749566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(osm->is_local[i], &idx_local)); 3759566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(osm->y[i], isll, osm->ly, isll_local, &osm->lprolongation[i])); 376d8b3b5e3Seaulisa 3779566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isll)); 3789566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isll_local)); 379d8b3b5e3Seaulisa } 38080ec0b7dSPatrick Sanan } 3819566063dSJacob Faibussowitsch PetscCall(VecDestroy(&vec)); 38280ec0b7dSPatrick Sanan } 38380ec0b7dSPatrick Sanan 384fb745f2cSMatthew G. Knepley if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE) { 385235cc4ceSMatthew G. Knepley IS *cis; 386235cc4ceSMatthew G. Knepley PetscInt c; 387235cc4ceSMatthew G. Knepley 3889566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(osm->n_local_true, &cis)); 389235cc4ceSMatthew G. Knepley for (c = 0; c < osm->n_local_true; ++c) cis[c] = osm->lis; 3909566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(pc->pmat, osm->n_local_true, osm->is, cis, scall, &osm->lmats)); 3919566063dSJacob Faibussowitsch PetscCall(PetscFree(cis)); 392fb745f2cSMatthew G. Knepley } 3934b9ad928SBarry Smith 3944b9ad928SBarry Smith /* Return control to the user so that the submatrices can be modified (e.g., to apply 3954b9ad928SBarry Smith different boundary conditions for the submatrices than for the global problem) */ 3969566063dSJacob Faibussowitsch PetscCall(PCModifySubMatrices(pc, osm->n_local_true, osm->is, osm->is, osm->pmat, pc->modifysubmatricesP)); 3974b9ad928SBarry Smith 39892bb6962SLisandro Dalcin /* 39992bb6962SLisandro Dalcin Loop over subdomains putting them into local ksp 40092bb6962SLisandro Dalcin */ 4019566063dSJacob Faibussowitsch PetscCall(KSPGetOptionsPrefix(osm->ksp[0], &prefix)); 40292bb6962SLisandro Dalcin for (i = 0; i < osm->n_local_true; i++) { 4039566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(osm->ksp[i], osm->pmat[i], osm->pmat[i])); 4049566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(osm->pmat[i], prefix)); 40548a46eb9SPierre Jolivet if (!pc->setupcalled) PetscCall(KSPSetFromOptions(osm->ksp[i])); 40692bb6962SLisandro Dalcin } 4073ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4084b9ad928SBarry Smith } 4094b9ad928SBarry Smith 410d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetUpOnBlocks_ASM(PC pc) 411d71ae5a4SJacob Faibussowitsch { 4124b9ad928SBarry Smith PC_ASM *osm = (PC_ASM *)pc->data; 41313f74950SBarry Smith PetscInt i; 414539c167fSBarry Smith KSPConvergedReason reason; 4154b9ad928SBarry Smith 4164b9ad928SBarry Smith PetscFunctionBegin; 4174b9ad928SBarry Smith for (i = 0; i < osm->n_local_true; i++) { 4189566063dSJacob Faibussowitsch PetscCall(KSPSetUp(osm->ksp[i])); 4199566063dSJacob Faibussowitsch PetscCall(KSPGetConvergedReason(osm->ksp[i], &reason)); 420ad540459SPierre Jolivet if (reason == KSP_DIVERGED_PC_FAILED) pc->failedreason = PC_SUBPC_ERROR; 4214b9ad928SBarry Smith } 4223ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4234b9ad928SBarry Smith } 4244b9ad928SBarry Smith 425d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApply_ASM(PC pc, Vec x, Vec y) 426d71ae5a4SJacob Faibussowitsch { 4274b9ad928SBarry Smith PC_ASM *osm = (PC_ASM *)pc->data; 4281dd8081eSeaulisa PetscInt i, n_local_true = osm->n_local_true; 4294b9ad928SBarry Smith ScatterMode forward = SCATTER_FORWARD, reverse = SCATTER_REVERSE; 4304b9ad928SBarry Smith 4314b9ad928SBarry Smith PetscFunctionBegin; 4324b9ad928SBarry Smith /* 43348e38a8aSPierre Jolivet support for limiting the restriction or interpolation to only local 4344b9ad928SBarry Smith subdomain values (leaving the other values 0). 4354b9ad928SBarry Smith */ 4364b9ad928SBarry Smith if (!(osm->type & PC_ASM_RESTRICT)) { 4374b9ad928SBarry Smith forward = SCATTER_FORWARD_LOCAL; 4384b9ad928SBarry Smith /* have to zero the work RHS since scatter may leave some slots empty */ 4399566063dSJacob Faibussowitsch PetscCall(VecSet(osm->lx, 0.0)); 4404b9ad928SBarry Smith } 441ad540459SPierre Jolivet if (!(osm->type & PC_ASM_INTERPOLATE)) reverse = SCATTER_REVERSE_LOCAL; 4424b9ad928SBarry Smith 4430fdf79fbSJacob 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]); 444b0de9f23SBarry Smith /* zero the global and the local solutions */ 4459566063dSJacob Faibussowitsch PetscCall(VecSet(y, 0.0)); 4469566063dSJacob Faibussowitsch PetscCall(VecSet(osm->ly, 0.0)); 447347574c9Seaulisa 44848e38a8aSPierre Jolivet /* copy the global RHS to local RHS including the ghost nodes */ 4499566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->restriction, x, osm->lx, INSERT_VALUES, forward)); 4509566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->restriction, x, osm->lx, INSERT_VALUES, forward)); 451347574c9Seaulisa 45248e38a8aSPierre Jolivet /* restrict local RHS to the overlapping 0-block RHS */ 4539566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward)); 4549566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward)); 455d8b3b5e3Seaulisa 45612cd4985SMatthew G. Knepley /* do the local solves */ 45712cd4985SMatthew G. Knepley for (i = 0; i < n_local_true; ++i) { 458b0de9f23SBarry Smith /* solve the overlapping i-block */ 4599566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PC_ApplyOnBlocks, osm->ksp[i], osm->x[i], osm->y[i], 0)); 4609566063dSJacob Faibussowitsch PetscCall(KSPSolve(osm->ksp[i], osm->x[i], osm->y[i])); 4619566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(osm->ksp[i], pc, osm->y[i])); 4629566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PC_ApplyOnBlocks, osm->ksp[i], osm->x[i], osm->y[i], 0)); 463d8b3b5e3Seaulisa 464677c7726SPierre Jolivet if (osm->lprolongation && osm->type != PC_ASM_INTERPOLATE) { /* interpolate the non-overlapping i-block solution to the local solution (only for restrictive additive) */ 4659566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward)); 4669566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward)); 46748e38a8aSPierre Jolivet } else { /* interpolate the overlapping i-block solution to the local solution */ 4689566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse)); 4699566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse)); 470d8b3b5e3Seaulisa } 471347574c9Seaulisa 472347574c9Seaulisa if (i < n_local_true - 1) { 47348e38a8aSPierre Jolivet /* restrict local RHS to the overlapping (i+1)-block RHS */ 4749566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->lrestriction[i + 1], osm->lx, osm->x[i + 1], INSERT_VALUES, forward)); 4759566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->lrestriction[i + 1], osm->lx, osm->x[i + 1], INSERT_VALUES, forward)); 476347574c9Seaulisa 477347574c9Seaulisa if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE) { 4788966356dSPierre Jolivet /* update the overlapping (i+1)-block RHS using the current local solution */ 4799566063dSJacob Faibussowitsch PetscCall(MatMult(osm->lmats[i + 1], osm->ly, osm->y[i + 1])); 4809566063dSJacob Faibussowitsch PetscCall(VecAXPBY(osm->x[i + 1], -1., 1., osm->y[i + 1])); 4817c3d802fSMatthew G. Knepley } 48212cd4985SMatthew G. Knepley } 48312cd4985SMatthew G. Knepley } 48448e38a8aSPierre Jolivet /* add the local solution to the global solution including the ghost nodes */ 4859566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->restriction, osm->ly, y, ADD_VALUES, reverse)); 4869566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->restriction, osm->ly, y, ADD_VALUES, reverse)); 4873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4884b9ad928SBarry Smith } 4894b9ad928SBarry Smith 490d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCMatApply_ASM(PC pc, Mat X, Mat Y) 491d71ae5a4SJacob Faibussowitsch { 49248e38a8aSPierre Jolivet PC_ASM *osm = (PC_ASM *)pc->data; 49348e38a8aSPierre Jolivet Mat Z, W; 49448e38a8aSPierre Jolivet Vec x; 49548e38a8aSPierre Jolivet PetscInt i, m, N; 49648e38a8aSPierre Jolivet ScatterMode forward = SCATTER_FORWARD, reverse = SCATTER_REVERSE; 49748e38a8aSPierre Jolivet 49848e38a8aSPierre Jolivet PetscFunctionBegin; 4997827d75bSBarry Smith PetscCheck(osm->n_local_true <= 1, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Not yet implemented"); 50048e38a8aSPierre Jolivet /* 50148e38a8aSPierre Jolivet support for limiting the restriction or interpolation to only local 50248e38a8aSPierre Jolivet subdomain values (leaving the other values 0). 50348e38a8aSPierre Jolivet */ 50448e38a8aSPierre Jolivet if (!(osm->type & PC_ASM_RESTRICT)) { 50548e38a8aSPierre Jolivet forward = SCATTER_FORWARD_LOCAL; 50648e38a8aSPierre Jolivet /* have to zero the work RHS since scatter may leave some slots empty */ 5079566063dSJacob Faibussowitsch PetscCall(VecSet(osm->lx, 0.0)); 50848e38a8aSPierre Jolivet } 509ad540459SPierre Jolivet if (!(osm->type & PC_ASM_INTERPOLATE)) reverse = SCATTER_REVERSE_LOCAL; 5109566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(osm->x[0], &m)); 5119566063dSJacob Faibussowitsch PetscCall(MatGetSize(X, NULL, &N)); 5129566063dSJacob Faibussowitsch PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, m, N, NULL, &Z)); 5130fdf79fbSJacob Faibussowitsch 5140fdf79fbSJacob 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]); 51548e38a8aSPierre Jolivet /* zero the global and the local solutions */ 5169566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(Y)); 5179566063dSJacob Faibussowitsch PetscCall(VecSet(osm->ly, 0.0)); 51848e38a8aSPierre Jolivet 51948e38a8aSPierre Jolivet for (i = 0; i < N; ++i) { 5209566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumnVecRead(X, i, &x)); 52148e38a8aSPierre Jolivet /* copy the global RHS to local RHS including the ghost nodes */ 5229566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->restriction, x, osm->lx, INSERT_VALUES, forward)); 5239566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->restriction, x, osm->lx, INSERT_VALUES, forward)); 5249566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumnVecRead(X, i, &x)); 52548e38a8aSPierre Jolivet 5269566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumnVecWrite(Z, i, &x)); 52748e38a8aSPierre Jolivet /* restrict local RHS to the overlapping 0-block RHS */ 5289566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->lrestriction[0], osm->lx, x, INSERT_VALUES, forward)); 5299566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->lrestriction[0], osm->lx, x, INSERT_VALUES, forward)); 5309566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumnVecWrite(Z, i, &x)); 53148e38a8aSPierre Jolivet } 5329566063dSJacob Faibussowitsch PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, m, N, NULL, &W)); 53348e38a8aSPierre Jolivet /* solve the overlapping 0-block */ 5349566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PC_ApplyOnBlocks, osm->ksp[0], Z, W, 0)); 5359566063dSJacob Faibussowitsch PetscCall(KSPMatSolve(osm->ksp[0], Z, W)); 5369566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(osm->ksp[0], pc, NULL)); 5379566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PC_ApplyOnBlocks, osm->ksp[0], Z, W, 0)); 5389566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Z)); 53948e38a8aSPierre Jolivet 54048e38a8aSPierre Jolivet for (i = 0; i < N; ++i) { 5419566063dSJacob Faibussowitsch PetscCall(VecSet(osm->ly, 0.0)); 5429566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumnVecRead(W, i, &x)); 543677c7726SPierre Jolivet if (osm->lprolongation && osm->type != PC_ASM_INTERPOLATE) { /* interpolate the non-overlapping 0-block solution to the local solution (only for restrictive additive) */ 5449566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->lprolongation[0], x, osm->ly, ADD_VALUES, forward)); 5459566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->lprolongation[0], x, osm->ly, ADD_VALUES, forward)); 54648e38a8aSPierre Jolivet } else { /* interpolate the overlapping 0-block solution to the local solution */ 5479566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->lrestriction[0], x, osm->ly, ADD_VALUES, reverse)); 5489566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->lrestriction[0], x, osm->ly, ADD_VALUES, reverse)); 54948e38a8aSPierre Jolivet } 5509566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumnVecRead(W, i, &x)); 55148e38a8aSPierre Jolivet 5529566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumnVecWrite(Y, i, &x)); 55348e38a8aSPierre Jolivet /* add the local solution to the global solution including the ghost nodes */ 5549566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->restriction, osm->ly, x, ADD_VALUES, reverse)); 5559566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->restriction, osm->ly, x, ADD_VALUES, reverse)); 5569566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumnVecWrite(Y, i, &x)); 55748e38a8aSPierre Jolivet } 5589566063dSJacob Faibussowitsch PetscCall(MatDestroy(&W)); 5593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 56048e38a8aSPierre Jolivet } 56148e38a8aSPierre Jolivet 562d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCApplyTranspose_ASM(PC pc, Vec x, Vec y) 563d71ae5a4SJacob Faibussowitsch { 5644b9ad928SBarry Smith PC_ASM *osm = (PC_ASM *)pc->data; 5651dd8081eSeaulisa PetscInt i, n_local_true = osm->n_local_true; 5664b9ad928SBarry Smith ScatterMode forward = SCATTER_FORWARD, reverse = SCATTER_REVERSE; 5674b9ad928SBarry Smith 5684b9ad928SBarry Smith PetscFunctionBegin; 5694b9ad928SBarry Smith /* 5704b9ad928SBarry Smith Support for limiting the restriction or interpolation to only local 5714b9ad928SBarry Smith subdomain values (leaving the other values 0). 5724b9ad928SBarry Smith 5734b9ad928SBarry Smith Note: these are reversed from the PCApply_ASM() because we are applying the 5744b9ad928SBarry Smith transpose of the three terms 5754b9ad928SBarry Smith */ 576d8b3b5e3Seaulisa 5774b9ad928SBarry Smith if (!(osm->type & PC_ASM_INTERPOLATE)) { 5784b9ad928SBarry Smith forward = SCATTER_FORWARD_LOCAL; 5794b9ad928SBarry Smith /* have to zero the work RHS since scatter may leave some slots empty */ 5809566063dSJacob Faibussowitsch PetscCall(VecSet(osm->lx, 0.0)); 5814b9ad928SBarry Smith } 5822fa5cd67SKarl Rupp if (!(osm->type & PC_ASM_RESTRICT)) reverse = SCATTER_REVERSE_LOCAL; 5834b9ad928SBarry Smith 584b0de9f23SBarry Smith /* zero the global and the local solutions */ 5859566063dSJacob Faibussowitsch PetscCall(VecSet(y, 0.0)); 5869566063dSJacob Faibussowitsch PetscCall(VecSet(osm->ly, 0.0)); 587d8b3b5e3Seaulisa 588b0de9f23SBarry Smith /* Copy the global RHS to local RHS including the ghost nodes */ 5899566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->restriction, x, osm->lx, INSERT_VALUES, forward)); 5909566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->restriction, x, osm->lx, INSERT_VALUES, forward)); 591d8b3b5e3Seaulisa 592b0de9f23SBarry Smith /* Restrict local RHS to the overlapping 0-block RHS */ 5939566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward)); 5949566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward)); 595d8b3b5e3Seaulisa 5964b9ad928SBarry Smith /* do the local solves */ 597d8b3b5e3Seaulisa for (i = 0; i < n_local_true; ++i) { 598b0de9f23SBarry Smith /* solve the overlapping i-block */ 5999566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PC_ApplyOnBlocks, osm->ksp[i], osm->x[i], osm->y[i], 0)); 6009566063dSJacob Faibussowitsch PetscCall(KSPSolveTranspose(osm->ksp[i], osm->x[i], osm->y[i])); 6019566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(osm->ksp[i], pc, osm->y[i])); 6029566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PC_ApplyOnBlocks, osm->ksp[i], osm->x[i], osm->y[i], 0)); 603d8b3b5e3Seaulisa 604a681a0f1SPierre Jolivet if (osm->lprolongation && osm->type != PC_ASM_RESTRICT) { /* interpolate the non-overlapping i-block solution to the local solution */ 6059566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward)); 6069566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward)); 607910cf402Sprj- } else { /* interpolate the overlapping i-block solution to the local solution */ 6089566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse)); 6099566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse)); 6104b9ad928SBarry Smith } 611d8b3b5e3Seaulisa 612d8b3b5e3Seaulisa if (i < n_local_true - 1) { 613b0de9f23SBarry Smith /* Restrict local RHS to the overlapping (i+1)-block RHS */ 6149566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->lrestriction[i + 1], osm->lx, osm->x[i + 1], INSERT_VALUES, forward)); 6159566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->lrestriction[i + 1], osm->lx, osm->x[i + 1], INSERT_VALUES, forward)); 6164b9ad928SBarry Smith } 6174b9ad928SBarry Smith } 618b0de9f23SBarry Smith /* Add the local solution to the global solution including the ghost nodes */ 6199566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->restriction, osm->ly, y, ADD_VALUES, reverse)); 6209566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->restriction, osm->ly, y, ADD_VALUES, reverse)); 6213ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6224b9ad928SBarry Smith } 6234b9ad928SBarry Smith 624d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCReset_ASM(PC pc) 625d71ae5a4SJacob Faibussowitsch { 6264b9ad928SBarry Smith PC_ASM *osm = (PC_ASM *)pc->data; 62713f74950SBarry Smith PetscInt i; 6284b9ad928SBarry Smith 6294b9ad928SBarry Smith PetscFunctionBegin; 63092bb6962SLisandro Dalcin if (osm->ksp) { 63148a46eb9SPierre Jolivet for (i = 0; i < osm->n_local_true; i++) PetscCall(KSPReset(osm->ksp[i])); 63292bb6962SLisandro Dalcin } 633e09e08ccSBarry Smith if (osm->pmat) { 63448a46eb9SPierre Jolivet if (osm->n_local_true > 0) PetscCall(MatDestroySubMatrices(osm->n_local_true, &osm->pmat)); 63592bb6962SLisandro Dalcin } 6361dd8081eSeaulisa if (osm->lrestriction) { 6379566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&osm->restriction)); 6381dd8081eSeaulisa for (i = 0; i < osm->n_local_true; i++) { 6399566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&osm->lrestriction[i])); 6409566063dSJacob Faibussowitsch if (osm->lprolongation) PetscCall(VecScatterDestroy(&osm->lprolongation[i])); 6419566063dSJacob Faibussowitsch PetscCall(VecDestroy(&osm->x[i])); 6429566063dSJacob Faibussowitsch PetscCall(VecDestroy(&osm->y[i])); 6434b9ad928SBarry Smith } 6449566063dSJacob Faibussowitsch PetscCall(PetscFree(osm->lrestriction)); 6459566063dSJacob Faibussowitsch if (osm->lprolongation) PetscCall(PetscFree(osm->lprolongation)); 6469566063dSJacob Faibussowitsch PetscCall(PetscFree(osm->x)); 6479566063dSJacob Faibussowitsch PetscCall(PetscFree(osm->y)); 64892bb6962SLisandro Dalcin } 6499566063dSJacob Faibussowitsch PetscCall(PCASMDestroySubdomains(osm->n_local_true, osm->is, osm->is_local)); 6509566063dSJacob Faibussowitsch PetscCall(ISDestroy(&osm->lis)); 6519566063dSJacob Faibussowitsch PetscCall(VecDestroy(&osm->lx)); 6529566063dSJacob Faibussowitsch PetscCall(VecDestroy(&osm->ly)); 65348a46eb9SPierre Jolivet if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE) PetscCall(MatDestroyMatrices(osm->n_local_true, &osm->lmats)); 6542fa5cd67SKarl Rupp 6559566063dSJacob Faibussowitsch PetscCall(PetscFree(osm->sub_mat_type)); 65680ec0b7dSPatrick Sanan 6570a545947SLisandro Dalcin osm->is = NULL; 6580a545947SLisandro Dalcin osm->is_local = NULL; 6593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 660e91c6855SBarry Smith } 661e91c6855SBarry Smith 662d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCDestroy_ASM(PC pc) 663d71ae5a4SJacob Faibussowitsch { 664e91c6855SBarry Smith PC_ASM *osm = (PC_ASM *)pc->data; 665e91c6855SBarry Smith PetscInt i; 666e91c6855SBarry Smith 667e91c6855SBarry Smith PetscFunctionBegin; 6689566063dSJacob Faibussowitsch PetscCall(PCReset_ASM(pc)); 669e91c6855SBarry Smith if (osm->ksp) { 67048a46eb9SPierre Jolivet for (i = 0; i < osm->n_local_true; i++) PetscCall(KSPDestroy(&osm->ksp[i])); 6719566063dSJacob Faibussowitsch PetscCall(PetscFree(osm->ksp)); 672e91c6855SBarry Smith } 6739566063dSJacob Faibussowitsch PetscCall(PetscFree(pc->data)); 67496322394SPierre Jolivet 6759566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetLocalSubdomains_C", NULL)); 6769566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetTotalSubdomains_C", NULL)); 6779566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetOverlap_C", NULL)); 6789566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetType_C", NULL)); 6799566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMGetType_C", NULL)); 6809566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetLocalType_C", NULL)); 6819566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMGetLocalType_C", NULL)); 6829566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetSortIndices_C", NULL)); 6839566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMGetSubKSP_C", NULL)); 6849566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMGetSubMatType_C", NULL)); 6859566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetSubMatType_C", NULL)); 6863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 6874b9ad928SBarry Smith } 6884b9ad928SBarry Smith 689d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCSetFromOptions_ASM(PC pc, PetscOptionItems *PetscOptionsObject) 690d71ae5a4SJacob Faibussowitsch { 6914b9ad928SBarry Smith PC_ASM *osm = (PC_ASM *)pc->data; 6929dcbbd2bSBarry Smith PetscInt blocks, ovl; 69357501b6eSBarry Smith PetscBool flg; 69492bb6962SLisandro Dalcin PCASMType asmtype; 69512cd4985SMatthew G. Knepley PCCompositeType loctype; 69680ec0b7dSPatrick Sanan char sub_mat_type[256]; 6974b9ad928SBarry Smith 6984b9ad928SBarry Smith PetscFunctionBegin; 699d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "Additive Schwarz options"); 7009566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_asm_dm_subdomains", "Use DMCreateDomainDecomposition() to define subdomains", "PCASMSetDMSubdomains", osm->dm_subdomains, &osm->dm_subdomains, &flg)); 7019566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_asm_blocks", "Number of subdomains", "PCASMSetTotalSubdomains", osm->n, &blocks, &flg)); 70265db9045SDmitry Karpeev if (flg) { 7039566063dSJacob Faibussowitsch PetscCall(PCASMSetTotalSubdomains(pc, blocks, NULL, NULL)); 704d709ea83SDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 70565db9045SDmitry Karpeev } 7069566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_asm_local_blocks", "Number of local subdomains", "PCASMSetLocalSubdomains", osm->n_local_true, &blocks, &flg)); 707342c94f9SMatthew G. Knepley if (flg) { 7089566063dSJacob Faibussowitsch PetscCall(PCASMSetLocalSubdomains(pc, blocks, NULL, NULL)); 709342c94f9SMatthew G. Knepley osm->dm_subdomains = PETSC_FALSE; 710342c94f9SMatthew G. Knepley } 7119566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_asm_overlap", "Number of grid points overlap", "PCASMSetOverlap", osm->overlap, &ovl, &flg)); 71265db9045SDmitry Karpeev if (flg) { 7139566063dSJacob Faibussowitsch PetscCall(PCASMSetOverlap(pc, ovl)); 714d709ea83SDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 71565db9045SDmitry Karpeev } 71690d69ab7SBarry Smith flg = PETSC_FALSE; 7179566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-pc_asm_type", "Type of restriction/extension", "PCASMSetType", PCASMTypes, (PetscEnum)osm->type, (PetscEnum *)&asmtype, &flg)); 7189566063dSJacob Faibussowitsch if (flg) PetscCall(PCASMSetType(pc, asmtype)); 71912cd4985SMatthew G. Knepley flg = PETSC_FALSE; 7209566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-pc_asm_local_type", "Type of local solver composition", "PCASMSetLocalType", PCCompositeTypes, (PetscEnum)osm->loctype, (PetscEnum *)&loctype, &flg)); 7219566063dSJacob Faibussowitsch if (flg) PetscCall(PCASMSetLocalType(pc, loctype)); 7229566063dSJacob Faibussowitsch PetscCall(PetscOptionsFList("-pc_asm_sub_mat_type", "Subsolve Matrix Type", "PCASMSetSubMatType", MatList, NULL, sub_mat_type, 256, &flg)); 7231baa6e33SBarry Smith if (flg) PetscCall(PCASMSetSubMatType(pc, sub_mat_type)); 724d0609cedSBarry Smith PetscOptionsHeadEnd(); 7253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7264b9ad928SBarry Smith } 7274b9ad928SBarry Smith 728d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCASMSetLocalSubdomains_ASM(PC pc, PetscInt n, IS is[], IS is_local[]) 729d71ae5a4SJacob Faibussowitsch { 7304b9ad928SBarry Smith PC_ASM *osm = (PC_ASM *)pc->data; 73192bb6962SLisandro Dalcin PetscInt i; 7324b9ad928SBarry Smith 7334b9ad928SBarry Smith PetscFunctionBegin; 73463a3b9bcSJacob Faibussowitsch PetscCheck(n >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Each process must have 1 or more blocks, n = %" PetscInt_FMT, n); 7357827d75bSBarry Smith PetscCheck(!pc->setupcalled || (n == osm->n_local_true && !is), PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "PCASMSetLocalSubdomains() should be called before calling PCSetUp()."); 736e7e72b3dSBarry Smith 7374b9ad928SBarry Smith if (!pc->setupcalled) { 73892bb6962SLisandro Dalcin if (is) { 7399566063dSJacob Faibussowitsch for (i = 0; i < n; i++) PetscCall(PetscObjectReference((PetscObject)is[i])); 74092bb6962SLisandro Dalcin } 741832fc9a5SMatthew Knepley if (is_local) { 7429566063dSJacob Faibussowitsch for (i = 0; i < n; i++) PetscCall(PetscObjectReference((PetscObject)is_local[i])); 743832fc9a5SMatthew Knepley } 7449566063dSJacob Faibussowitsch PetscCall(PCASMDestroySubdomains(osm->n_local_true, osm->is, osm->is_local)); 7452fa5cd67SKarl Rupp 74607517c86SMark Adams if (osm->ksp && osm->n_local_true != n) { 74707517c86SMark Adams for (i = 0; i < osm->n_local_true; i++) PetscCall(KSPDestroy(&osm->ksp[i])); 74807517c86SMark Adams PetscCall(PetscFree(osm->ksp)); 74907517c86SMark Adams } 75007517c86SMark Adams 7514b9ad928SBarry Smith osm->n_local_true = n; 7520a545947SLisandro Dalcin osm->is = NULL; 7530a545947SLisandro Dalcin osm->is_local = NULL; 75492bb6962SLisandro Dalcin if (is) { 7559566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &osm->is)); 7562fa5cd67SKarl Rupp for (i = 0; i < n; i++) osm->is[i] = is[i]; 7573d03bb48SJed Brown /* Flag indicating that the user has set overlapping subdomains so PCASM should not increase their size. */ 7583d03bb48SJed Brown osm->overlap = -1; 75992bb6962SLisandro Dalcin } 7602b691e39SMatthew Knepley if (is_local) { 7619566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &osm->is_local)); 7622fa5cd67SKarl Rupp for (i = 0; i < n; i++) osm->is_local[i] = is_local[i]; 763a35b7b57SDmitry Karpeev if (!is) { 7649566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(osm->n_local_true, &osm->is)); 765a35b7b57SDmitry Karpeev for (i = 0; i < osm->n_local_true; i++) { 766a35b7b57SDmitry Karpeev if (osm->overlap > 0) { /* With positive overlap, osm->is[i] will be modified */ 7679566063dSJacob Faibussowitsch PetscCall(ISDuplicate(osm->is_local[i], &osm->is[i])); 7689566063dSJacob Faibussowitsch PetscCall(ISCopy(osm->is_local[i], osm->is[i])); 769a35b7b57SDmitry Karpeev } else { 7709566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)osm->is_local[i])); 771a35b7b57SDmitry Karpeev osm->is[i] = osm->is_local[i]; 772a35b7b57SDmitry Karpeev } 773a35b7b57SDmitry Karpeev } 774a35b7b57SDmitry Karpeev } 7752b691e39SMatthew Knepley } 7764b9ad928SBarry Smith } 7773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 7784b9ad928SBarry Smith } 7794b9ad928SBarry Smith 780d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCASMSetTotalSubdomains_ASM(PC pc, PetscInt N, IS *is, IS *is_local) 781d71ae5a4SJacob Faibussowitsch { 7824b9ad928SBarry Smith PC_ASM *osm = (PC_ASM *)pc->data; 78313f74950SBarry Smith PetscMPIInt rank, size; 78478904715SLisandro Dalcin PetscInt n; 7854b9ad928SBarry Smith 7864b9ad928SBarry Smith PetscFunctionBegin; 78763a3b9bcSJacob Faibussowitsch PetscCheck(N >= 1, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Number of total blocks must be > 0, N = %" PetscInt_FMT, N); 7887827d75bSBarry Smith PetscCheck(!is && !is_local, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Use PCASMSetLocalSubdomains() to set specific index sets\n\they cannot be set globally yet."); 7894b9ad928SBarry Smith 7904b9ad928SBarry Smith /* 791880770e9SJed Brown Split the subdomains equally among all processors 7924b9ad928SBarry Smith */ 7939566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank)); 7949566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size)); 7954b9ad928SBarry Smith n = N / size + ((N % size) > rank); 79663a3b9bcSJacob 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); 7977827d75bSBarry Smith PetscCheck(!pc->setupcalled || n == osm->n_local_true, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "PCASMSetTotalSubdomains() should be called before PCSetUp()."); 7984b9ad928SBarry Smith if (!pc->setupcalled) { 7999566063dSJacob Faibussowitsch PetscCall(PCASMDestroySubdomains(osm->n_local_true, osm->is, osm->is_local)); 8002fa5cd67SKarl Rupp 8014b9ad928SBarry Smith osm->n_local_true = n; 8020a545947SLisandro Dalcin osm->is = NULL; 8030a545947SLisandro Dalcin osm->is_local = NULL; 8044b9ad928SBarry Smith } 8053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8064b9ad928SBarry Smith } 8074b9ad928SBarry Smith 808d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCASMSetOverlap_ASM(PC pc, PetscInt ovl) 809d71ae5a4SJacob Faibussowitsch { 81092bb6962SLisandro Dalcin PC_ASM *osm = (PC_ASM *)pc->data; 8114b9ad928SBarry Smith 8124b9ad928SBarry Smith PetscFunctionBegin; 8137827d75bSBarry Smith PetscCheck(ovl >= 0, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Negative overlap value requested"); 8147827d75bSBarry Smith PetscCheck(!pc->setupcalled || ovl == osm->overlap, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "PCASMSetOverlap() should be called before PCSetUp()."); 8152fa5cd67SKarl Rupp if (!pc->setupcalled) osm->overlap = ovl; 8163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8174b9ad928SBarry Smith } 8184b9ad928SBarry Smith 819d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCASMSetType_ASM(PC pc, PCASMType type) 820d71ae5a4SJacob Faibussowitsch { 82192bb6962SLisandro Dalcin PC_ASM *osm = (PC_ASM *)pc->data; 8224b9ad928SBarry Smith 8234b9ad928SBarry Smith PetscFunctionBegin; 8244b9ad928SBarry Smith osm->type = type; 825bf108f30SBarry Smith osm->type_set = PETSC_TRUE; 8263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8274b9ad928SBarry Smith } 8284b9ad928SBarry Smith 829d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCASMGetType_ASM(PC pc, PCASMType *type) 830d71ae5a4SJacob Faibussowitsch { 831c60c7ad4SBarry Smith PC_ASM *osm = (PC_ASM *)pc->data; 832c60c7ad4SBarry Smith 833c60c7ad4SBarry Smith PetscFunctionBegin; 834c60c7ad4SBarry Smith *type = osm->type; 8353ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 836c60c7ad4SBarry Smith } 837c60c7ad4SBarry Smith 838d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCASMSetLocalType_ASM(PC pc, PCCompositeType type) 839d71ae5a4SJacob Faibussowitsch { 84012cd4985SMatthew G. Knepley PC_ASM *osm = (PC_ASM *)pc->data; 84112cd4985SMatthew G. Knepley 84212cd4985SMatthew G. Knepley PetscFunctionBegin; 8437827d75bSBarry 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"); 84412cd4985SMatthew G. Knepley osm->loctype = type; 8453ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 84612cd4985SMatthew G. Knepley } 84712cd4985SMatthew G. Knepley 848d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCASMGetLocalType_ASM(PC pc, PCCompositeType *type) 849d71ae5a4SJacob Faibussowitsch { 85012cd4985SMatthew G. Knepley PC_ASM *osm = (PC_ASM *)pc->data; 85112cd4985SMatthew G. Knepley 85212cd4985SMatthew G. Knepley PetscFunctionBegin; 85312cd4985SMatthew G. Knepley *type = osm->loctype; 8543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 85512cd4985SMatthew G. Knepley } 85612cd4985SMatthew G. Knepley 857d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCASMSetSortIndices_ASM(PC pc, PetscBool doSort) 858d71ae5a4SJacob Faibussowitsch { 8596ed231c7SMatthew Knepley PC_ASM *osm = (PC_ASM *)pc->data; 8606ed231c7SMatthew Knepley 8616ed231c7SMatthew Knepley PetscFunctionBegin; 8626ed231c7SMatthew Knepley osm->sort_indices = doSort; 8633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8646ed231c7SMatthew Knepley } 8656ed231c7SMatthew Knepley 866d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCASMGetSubKSP_ASM(PC pc, PetscInt *n_local, PetscInt *first_local, KSP **ksp) 867d71ae5a4SJacob Faibussowitsch { 86892bb6962SLisandro Dalcin PC_ASM *osm = (PC_ASM *)pc->data; 8694b9ad928SBarry Smith 8704b9ad928SBarry Smith PetscFunctionBegin; 8717827d75bSBarry 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"); 8724b9ad928SBarry Smith 8732fa5cd67SKarl Rupp if (n_local) *n_local = osm->n_local_true; 87492bb6962SLisandro Dalcin if (first_local) { 8759566063dSJacob Faibussowitsch PetscCallMPI(MPI_Scan(&osm->n_local_true, first_local, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)pc))); 87692bb6962SLisandro Dalcin *first_local -= osm->n_local_true; 87792bb6962SLisandro Dalcin } 878ed155784SPierre Jolivet if (ksp) *ksp = osm->ksp; 8793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 8804b9ad928SBarry Smith } 8814b9ad928SBarry Smith 882d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCASMGetSubMatType_ASM(PC pc, MatType *sub_mat_type) 883d71ae5a4SJacob Faibussowitsch { 88480ec0b7dSPatrick Sanan PC_ASM *osm = (PC_ASM *)pc->data; 88580ec0b7dSPatrick Sanan 88680ec0b7dSPatrick Sanan PetscFunctionBegin; 88780ec0b7dSPatrick Sanan PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 8884f572ea9SToby Isaac PetscAssertPointer(sub_mat_type, 2); 88980ec0b7dSPatrick Sanan *sub_mat_type = osm->sub_mat_type; 8903ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 89180ec0b7dSPatrick Sanan } 89280ec0b7dSPatrick Sanan 893d71ae5a4SJacob Faibussowitsch static PetscErrorCode PCASMSetSubMatType_ASM(PC pc, MatType sub_mat_type) 894d71ae5a4SJacob Faibussowitsch { 89580ec0b7dSPatrick Sanan PC_ASM *osm = (PC_ASM *)pc->data; 89680ec0b7dSPatrick Sanan 89780ec0b7dSPatrick Sanan PetscFunctionBegin; 89880ec0b7dSPatrick Sanan PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 8999566063dSJacob Faibussowitsch PetscCall(PetscFree(osm->sub_mat_type)); 9009566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(sub_mat_type, (char **)&osm->sub_mat_type)); 9013ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 90280ec0b7dSPatrick Sanan } 90380ec0b7dSPatrick Sanan 9044b9ad928SBarry Smith /*@C 905f1580f4eSBarry Smith PCASMSetLocalSubdomains - Sets the local subdomains (for this processor only) for the additive Schwarz preconditioner `PCASM`. 9064b9ad928SBarry Smith 907c3339decSBarry Smith Collective 9084b9ad928SBarry Smith 9094b9ad928SBarry Smith Input Parameters: 9104b9ad928SBarry Smith + pc - the preconditioner context 9114b9ad928SBarry Smith . n - the number of subdomains for this processor (default value = 1) 9128c03b21aSDmitry Karpeev . is - the index set that defines the subdomains for this processor 91320f4b53cSBarry Smith (or `NULL` for PETSc to determine subdomains) 914f1ee410cSBarry 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 91520f4b53cSBarry Smith (or `NULL` to not provide these) 9164b9ad928SBarry Smith 917342c94f9SMatthew G. Knepley Options Database Key: 91820f4b53cSBarry Smith . -pc_asm_local_blocks <blks> - Sets number of local blocks 91920f4b53cSBarry Smith 92020f4b53cSBarry Smith Level: advanced 921342c94f9SMatthew G. Knepley 9224b9ad928SBarry Smith Notes: 923f1580f4eSBarry Smith The `IS` numbering is in the parallel, global numbering of the vector for both is and is_local 9244b9ad928SBarry Smith 925f1580f4eSBarry Smith By default the `PCASM` preconditioner uses 1 block per processor. 9264b9ad928SBarry Smith 927f1580f4eSBarry Smith Use `PCASMSetTotalSubdomains()` to set the subdomains for all processors. 9284b9ad928SBarry Smith 929f1580f4eSBarry Smith If is_local is provided and `PCASMType` is `PC_ASM_RESTRICT` then the solution only over the is_local region is interpolated 930f1ee410cSBarry Smith back to form the global solution (this is the standard restricted additive Schwarz method) 931f1ee410cSBarry Smith 932f1580f4eSBarry Smith If the is_local is provided and `PCASMType` is `PC_ASM_INTERPOLATE` or `PC_ASM_NONE` then an error is generated since there is 933f1ee410cSBarry Smith no code to handle that case. 934f1ee410cSBarry Smith 935562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMSetTotalSubdomains()`, `PCASMSetOverlap()`, `PCASMGetSubKSP()`, 936f1580f4eSBarry Smith `PCASMCreateSubdomains2D()`, `PCASMGetLocalSubdomains()`, `PCASMType`, `PCASMSetType()`, `PCGASM` 9374b9ad928SBarry Smith @*/ 938d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMSetLocalSubdomains(PC pc, PetscInt n, IS is[], IS is_local[]) 939d71ae5a4SJacob Faibussowitsch { 9404b9ad928SBarry Smith PetscFunctionBegin; 9410700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 942cac4c232SBarry Smith PetscTryMethod(pc, "PCASMSetLocalSubdomains_C", (PC, PetscInt, IS[], IS[]), (pc, n, is, is_local)); 9433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 9444b9ad928SBarry Smith } 9454b9ad928SBarry Smith 9464b9ad928SBarry Smith /*@C 947feb221f8SDmitry Karpeev PCASMSetTotalSubdomains - Sets the subdomains for all processors for the 948f1580f4eSBarry Smith additive Schwarz preconditioner, `PCASM`. 9494b9ad928SBarry Smith 950c3339decSBarry Smith Collective, all MPI ranks must pass in the same array of `IS` 9514b9ad928SBarry Smith 9524b9ad928SBarry Smith Input Parameters: 9534b9ad928SBarry Smith + pc - the preconditioner context 954feb221f8SDmitry Karpeev . N - the number of subdomains for all processors 955feb221f8SDmitry Karpeev . is - the index sets that define the subdomains for all processors 95620f4b53cSBarry Smith (or `NULL` to ask PETSc to determine the subdomains) 9572b691e39SMatthew Knepley - is_local - the index sets that define the local part of the subdomains for this processor 95820f4b53cSBarry Smith (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: 966f1ee410cSBarry 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: 12000298fd71SBarry Smith + n_local - the number of blocks on this processor or NULL 12010298fd71SBarry Smith . first_local - the global number of the first block on this processor or NULL, 12020298fd71SBarry Smith all processors must request or all must pass NULL 1203f1580f4eSBarry Smith - ksp - the array of `KSP` contexts 12044b9ad928SBarry Smith 120520f4b53cSBarry Smith Level: advanced 120620f4b53cSBarry Smith 1207f1580f4eSBarry Smith Notes: 1208f1580f4eSBarry Smith After `PCASMGetSubKSP()` the array of `KSP`s is not to be freed. 12094b9ad928SBarry Smith 1210f1580f4eSBarry Smith You must call `KSPSetUp()` before calling `PCASMGetSubKSP()`. 12114b9ad928SBarry Smith 1212feefa0e1SJacob Faibussowitsch Fortran Notes: 1213f1580f4eSBarry 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. 1214d29017ddSJed Brown 1215562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMSetTotalSubdomains()`, `PCASMSetOverlap()`, 1216db781477SPatrick Sanan `PCASMCreateSubdomains2D()`, 12174b9ad928SBarry Smith @*/ 1218d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMGetSubKSP(PC pc, PetscInt *n_local, PetscInt *first_local, KSP *ksp[]) 1219d71ae5a4SJacob Faibussowitsch { 12204b9ad928SBarry Smith PetscFunctionBegin; 12210700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1222cac4c232SBarry Smith PetscUseMethod(pc, "PCASMGetSubKSP_C", (PC, PetscInt *, PetscInt *, KSP **), (pc, n_local, first_local, ksp)); 12233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 12244b9ad928SBarry Smith } 12254b9ad928SBarry Smith 12264b9ad928SBarry Smith /*MC 12273b09bd56SBarry Smith PCASM - Use the (restricted) additive Schwarz method, each block is (approximately) solved with 1228*1d27aa22SBarry Smith its own `KSP` object, {cite}`dryja1987additive` and {cite}`1sbg` 12294b9ad928SBarry Smith 12304b9ad928SBarry Smith Options Database Keys: 1231*1d27aa22SBarry Smith + -pc_asm_blocks <blks> - Sets total blocks. Defaults to one block per MPI process. 12324b9ad928SBarry Smith . -pc_asm_overlap <ovl> - Sets overlap 1233f1580f4eSBarry Smith . -pc_asm_type [basic,restrict,interpolate,none] - Sets `PCASMType`, default is restrict. See `PCASMSetType()` 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 1240f1580f4eSBarry Smith will get a different convergence rate due to the default option of -pc_asm_type restrict. Use 1241f1580f4eSBarry 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 1246f1580f4eSBarry Smith To set options on the solvers for each block append -sub_ to all the `KSP`, and `PC` 1247f1580f4eSBarry 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 1252562efe2eSBarry Smith .seealso: [](ch_ksp), `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCASMType`, `PCCompositeType`, 1253db781477SPatrick Sanan `PCBJACOBI`, `PCASMGetSubKSP()`, `PCASMSetLocalSubdomains()`, `PCASMType`, `PCASMGetType()`, `PCASMSetLocalType()`, `PCASMGetLocalType()` 1254db781477SPatrick Sanan `PCASMSetTotalSubdomains()`, `PCSetModifySubMatrices()`, `PCASMSetOverlap()`, `PCASMSetType()`, `PCCompositeType` 12554b9ad928SBarry Smith M*/ 12564b9ad928SBarry Smith 1257d71ae5a4SJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCCreate_ASM(PC pc) 1258d71ae5a4SJacob Faibussowitsch { 12594b9ad928SBarry Smith PC_ASM *osm; 12604b9ad928SBarry Smith 12614b9ad928SBarry Smith PetscFunctionBegin; 12624dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&osm)); 12632fa5cd67SKarl Rupp 12644b9ad928SBarry Smith osm->n = PETSC_DECIDE; 12654b9ad928SBarry Smith osm->n_local = 0; 12662b837212SDmitry Karpeev osm->n_local_true = PETSC_DECIDE; 12674b9ad928SBarry Smith osm->overlap = 1; 12680a545947SLisandro Dalcin osm->ksp = NULL; 12690a545947SLisandro Dalcin osm->restriction = NULL; 12700a545947SLisandro Dalcin osm->lprolongation = NULL; 12710a545947SLisandro Dalcin osm->lrestriction = NULL; 12720a545947SLisandro Dalcin osm->x = NULL; 12730a545947SLisandro Dalcin osm->y = NULL; 12740a545947SLisandro Dalcin osm->is = NULL; 12750a545947SLisandro Dalcin osm->is_local = NULL; 12760a545947SLisandro Dalcin osm->mat = NULL; 12770a545947SLisandro Dalcin osm->pmat = NULL; 12784b9ad928SBarry Smith osm->type = PC_ASM_RESTRICT; 127912cd4985SMatthew G. Knepley osm->loctype = PC_COMPOSITE_ADDITIVE; 12806ed231c7SMatthew Knepley osm->sort_indices = PETSC_TRUE; 1281d709ea83SDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 128280ec0b7dSPatrick Sanan osm->sub_mat_type = NULL; 12834b9ad928SBarry Smith 128492bb6962SLisandro Dalcin pc->data = (void *)osm; 12854b9ad928SBarry Smith pc->ops->apply = PCApply_ASM; 128648e38a8aSPierre Jolivet pc->ops->matapply = PCMatApply_ASM; 12874b9ad928SBarry Smith pc->ops->applytranspose = PCApplyTranspose_ASM; 12884b9ad928SBarry Smith pc->ops->setup = PCSetUp_ASM; 1289e91c6855SBarry Smith pc->ops->reset = PCReset_ASM; 12904b9ad928SBarry Smith pc->ops->destroy = PCDestroy_ASM; 12914b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_ASM; 12924b9ad928SBarry Smith pc->ops->setuponblocks = PCSetUpOnBlocks_ASM; 12934b9ad928SBarry Smith pc->ops->view = PCView_ASM; 12940a545947SLisandro Dalcin pc->ops->applyrichardson = NULL; 12954b9ad928SBarry Smith 12969566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetLocalSubdomains_C", PCASMSetLocalSubdomains_ASM)); 12979566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetTotalSubdomains_C", PCASMSetTotalSubdomains_ASM)); 12989566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetOverlap_C", PCASMSetOverlap_ASM)); 12999566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetType_C", PCASMSetType_ASM)); 13009566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMGetType_C", PCASMGetType_ASM)); 13019566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetLocalType_C", PCASMSetLocalType_ASM)); 13029566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMGetLocalType_C", PCASMGetLocalType_ASM)); 13039566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetSortIndices_C", PCASMSetSortIndices_ASM)); 13049566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMGetSubKSP_C", PCASMGetSubKSP_ASM)); 13059566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMGetSubMatType_C", PCASMGetSubMatType_ASM)); 13069566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetSubMatType_C", PCASMSetSubMatType_ASM)); 13073ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13084b9ad928SBarry Smith } 13094b9ad928SBarry Smith 131092bb6962SLisandro Dalcin /*@C 131192bb6962SLisandro Dalcin PCASMCreateSubdomains - Creates the index sets for the overlapping Schwarz 1312f1580f4eSBarry Smith preconditioner, `PCASM`, for any problem on a general grid. 131392bb6962SLisandro Dalcin 131492bb6962SLisandro Dalcin Collective 131592bb6962SLisandro Dalcin 131692bb6962SLisandro Dalcin Input Parameters: 131792bb6962SLisandro Dalcin + A - The global matrix operator 131892bb6962SLisandro Dalcin - n - the number of local blocks 131992bb6962SLisandro Dalcin 13202fe279fdSBarry Smith Output Parameter: 132192bb6962SLisandro Dalcin . outis - the array of index sets defining the subdomains 132292bb6962SLisandro Dalcin 132392bb6962SLisandro Dalcin Level: advanced 132492bb6962SLisandro Dalcin 1325f1580f4eSBarry Smith Note: 1326f1580f4eSBarry Smith This generates nonoverlapping subdomains; the `PCASM` will generate the overlap 1327f1580f4eSBarry Smith from these if you use `PCASMSetLocalSubdomains()` 13287d6bfa3bSBarry Smith 1329feefa0e1SJacob Faibussowitsch Fortran Notes: 1330f1580f4eSBarry Smith You must provide the array outis[] already allocated of length n. 13317d6bfa3bSBarry Smith 1332562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMSetLocalSubdomains()`, `PCASMDestroySubdomains()` 133392bb6962SLisandro Dalcin @*/ 1334d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMCreateSubdomains(Mat A, PetscInt n, IS *outis[]) 1335d71ae5a4SJacob Faibussowitsch { 133692bb6962SLisandro Dalcin MatPartitioning mpart; 133792bb6962SLisandro Dalcin const char *prefix; 133892bb6962SLisandro Dalcin PetscInt i, j, rstart, rend, bs; 1339976e8c5aSLisandro Dalcin PetscBool hasop, isbaij = PETSC_FALSE, foundpart = PETSC_FALSE; 13400298fd71SBarry Smith Mat Ad = NULL, adj; 134192bb6962SLisandro Dalcin IS ispart, isnumb, *is; 134292bb6962SLisandro Dalcin 134392bb6962SLisandro Dalcin PetscFunctionBegin; 13440700a824SBarry Smith PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 13454f572ea9SToby Isaac PetscAssertPointer(outis, 3); 134663a3b9bcSJacob Faibussowitsch PetscCheck(n >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "number of local blocks must be > 0, n = %" PetscInt_FMT, n); 134792bb6962SLisandro Dalcin 134892bb6962SLisandro Dalcin /* Get prefix, row distribution, and block size */ 13499566063dSJacob Faibussowitsch PetscCall(MatGetOptionsPrefix(A, &prefix)); 13509566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A, &rstart, &rend)); 13519566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(A, &bs)); 135263a3b9bcSJacob 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); 135365e19b50SBarry Smith 135492bb6962SLisandro Dalcin /* Get diagonal block from matrix if possible */ 13559566063dSJacob Faibussowitsch PetscCall(MatHasOperation(A, MATOP_GET_DIAGONAL_BLOCK, &hasop)); 135648a46eb9SPierre Jolivet if (hasop) PetscCall(MatGetDiagonalBlock(A, &Ad)); 135792bb6962SLisandro Dalcin if (Ad) { 13589566063dSJacob Faibussowitsch PetscCall(PetscObjectBaseTypeCompare((PetscObject)Ad, MATSEQBAIJ, &isbaij)); 13599566063dSJacob Faibussowitsch if (!isbaij) PetscCall(PetscObjectBaseTypeCompare((PetscObject)Ad, MATSEQSBAIJ, &isbaij)); 136092bb6962SLisandro Dalcin } 136192bb6962SLisandro Dalcin if (Ad && n > 1) { 1362ace3abfcSBarry Smith PetscBool match, done; 136392bb6962SLisandro Dalcin /* Try to setup a good matrix partitioning if available */ 13649566063dSJacob Faibussowitsch PetscCall(MatPartitioningCreate(PETSC_COMM_SELF, &mpart)); 13659566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)mpart, prefix)); 13669566063dSJacob Faibussowitsch PetscCall(MatPartitioningSetFromOptions(mpart)); 13679566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)mpart, MATPARTITIONINGCURRENT, &match)); 136848a46eb9SPierre Jolivet if (!match) PetscCall(PetscObjectTypeCompare((PetscObject)mpart, MATPARTITIONINGSQUARE, &match)); 136992bb6962SLisandro Dalcin if (!match) { /* assume a "good" partitioner is available */ 13701a83f524SJed Brown PetscInt na; 13711a83f524SJed Brown const PetscInt *ia, *ja; 13729566063dSJacob Faibussowitsch PetscCall(MatGetRowIJ(Ad, 0, PETSC_TRUE, isbaij, &na, &ia, &ja, &done)); 137392bb6962SLisandro Dalcin if (done) { 137492bb6962SLisandro Dalcin /* Build adjacency matrix by hand. Unfortunately a call to 137592bb6962SLisandro Dalcin MatConvert(Ad,MATMPIADJ,MAT_INITIAL_MATRIX,&adj) will 137692bb6962SLisandro Dalcin remove the block-aij structure and we cannot expect 137792bb6962SLisandro Dalcin MatPartitioning to split vertices as we need */ 13780a545947SLisandro Dalcin PetscInt i, j, len, nnz, cnt, *iia = NULL, *jja = NULL; 13791a83f524SJed Brown const PetscInt *row; 138092bb6962SLisandro Dalcin nnz = 0; 138192bb6962SLisandro Dalcin for (i = 0; i < na; i++) { /* count number of nonzeros */ 138292bb6962SLisandro Dalcin len = ia[i + 1] - ia[i]; 138392bb6962SLisandro Dalcin row = ja + ia[i]; 138492bb6962SLisandro Dalcin for (j = 0; j < len; j++) { 138592bb6962SLisandro Dalcin if (row[j] == i) { /* don't count diagonal */ 13869371c9d4SSatish Balay len--; 13879371c9d4SSatish Balay break; 138892bb6962SLisandro Dalcin } 138992bb6962SLisandro Dalcin } 139092bb6962SLisandro Dalcin nnz += len; 139192bb6962SLisandro Dalcin } 13929566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(na + 1, &iia)); 13939566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nnz, &jja)); 139492bb6962SLisandro Dalcin nnz = 0; 139592bb6962SLisandro Dalcin iia[0] = 0; 139692bb6962SLisandro Dalcin for (i = 0; i < na; i++) { /* fill adjacency */ 139792bb6962SLisandro Dalcin cnt = 0; 139892bb6962SLisandro Dalcin len = ia[i + 1] - ia[i]; 139992bb6962SLisandro Dalcin row = ja + ia[i]; 140092bb6962SLisandro Dalcin for (j = 0; j < len; j++) { 140192bb6962SLisandro Dalcin if (row[j] != i) { /* if not diagonal */ 140292bb6962SLisandro Dalcin jja[nnz + cnt++] = row[j]; 140392bb6962SLisandro Dalcin } 140492bb6962SLisandro Dalcin } 140592bb6962SLisandro Dalcin nnz += cnt; 140692bb6962SLisandro Dalcin iia[i + 1] = nnz; 140792bb6962SLisandro Dalcin } 140892bb6962SLisandro Dalcin /* Partitioning of the adjacency matrix */ 14099566063dSJacob Faibussowitsch PetscCall(MatCreateMPIAdj(PETSC_COMM_SELF, na, na, iia, jja, NULL, &adj)); 14109566063dSJacob Faibussowitsch PetscCall(MatPartitioningSetAdjacency(mpart, adj)); 14119566063dSJacob Faibussowitsch PetscCall(MatPartitioningSetNParts(mpart, n)); 14129566063dSJacob Faibussowitsch PetscCall(MatPartitioningApply(mpart, &ispart)); 14139566063dSJacob Faibussowitsch PetscCall(ISPartitioningToNumbering(ispart, &isnumb)); 14149566063dSJacob Faibussowitsch PetscCall(MatDestroy(&adj)); 141592bb6962SLisandro Dalcin foundpart = PETSC_TRUE; 141692bb6962SLisandro Dalcin } 14179566063dSJacob Faibussowitsch PetscCall(MatRestoreRowIJ(Ad, 0, PETSC_TRUE, isbaij, &na, &ia, &ja, &done)); 141892bb6962SLisandro Dalcin } 14199566063dSJacob Faibussowitsch PetscCall(MatPartitioningDestroy(&mpart)); 142092bb6962SLisandro Dalcin } 142192bb6962SLisandro Dalcin 14229566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &is)); 142392bb6962SLisandro Dalcin *outis = is; 142492bb6962SLisandro Dalcin 142592bb6962SLisandro Dalcin if (!foundpart) { 142692bb6962SLisandro Dalcin /* Partitioning by contiguous chunks of rows */ 142792bb6962SLisandro Dalcin 142892bb6962SLisandro Dalcin PetscInt mbs = (rend - rstart) / bs; 142992bb6962SLisandro Dalcin PetscInt start = rstart; 143092bb6962SLisandro Dalcin for (i = 0; i < n; i++) { 143192bb6962SLisandro Dalcin PetscInt count = (mbs / n + ((mbs % n) > i)) * bs; 14329566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, count, start, 1, &is[i])); 143392bb6962SLisandro Dalcin start += count; 143492bb6962SLisandro Dalcin } 143592bb6962SLisandro Dalcin 143692bb6962SLisandro Dalcin } else { 143792bb6962SLisandro Dalcin /* Partitioning by adjacency of diagonal block */ 143892bb6962SLisandro Dalcin 143992bb6962SLisandro Dalcin const PetscInt *numbering; 144092bb6962SLisandro Dalcin PetscInt *count, nidx, *indices, *newidx, start = 0; 144192bb6962SLisandro Dalcin /* Get node count in each partition */ 14429566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &count)); 14439566063dSJacob Faibussowitsch PetscCall(ISPartitioningCount(ispart, n, count)); 144492bb6962SLisandro Dalcin if (isbaij && bs > 1) { /* adjust for the block-aij case */ 144592bb6962SLisandro Dalcin for (i = 0; i < n; i++) count[i] *= bs; 144692bb6962SLisandro Dalcin } 144792bb6962SLisandro Dalcin /* Build indices from node numbering */ 14489566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(isnumb, &nidx)); 14499566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nidx, &indices)); 145092bb6962SLisandro Dalcin for (i = 0; i < nidx; i++) indices[i] = i; /* needs to be initialized */ 14519566063dSJacob Faibussowitsch PetscCall(ISGetIndices(isnumb, &numbering)); 14529566063dSJacob Faibussowitsch PetscCall(PetscSortIntWithPermutation(nidx, numbering, indices)); 14539566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(isnumb, &numbering)); 145492bb6962SLisandro Dalcin if (isbaij && bs > 1) { /* adjust for the block-aij case */ 14559566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nidx * bs, &newidx)); 14562fa5cd67SKarl Rupp for (i = 0; i < nidx; i++) { 14572fa5cd67SKarl Rupp for (j = 0; j < bs; j++) newidx[i * bs + j] = indices[i] * bs + j; 14582fa5cd67SKarl Rupp } 14599566063dSJacob Faibussowitsch PetscCall(PetscFree(indices)); 146092bb6962SLisandro Dalcin nidx *= bs; 146192bb6962SLisandro Dalcin indices = newidx; 146292bb6962SLisandro Dalcin } 146392bb6962SLisandro Dalcin /* Shift to get global indices */ 146492bb6962SLisandro Dalcin for (i = 0; i < nidx; i++) indices[i] += rstart; 146592bb6962SLisandro Dalcin 146692bb6962SLisandro Dalcin /* Build the index sets for each block */ 146792bb6962SLisandro Dalcin for (i = 0; i < n; i++) { 14689566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, count[i], &indices[start], PETSC_COPY_VALUES, &is[i])); 14699566063dSJacob Faibussowitsch PetscCall(ISSort(is[i])); 147092bb6962SLisandro Dalcin start += count[i]; 147192bb6962SLisandro Dalcin } 147292bb6962SLisandro Dalcin 14739566063dSJacob Faibussowitsch PetscCall(PetscFree(count)); 14749566063dSJacob Faibussowitsch PetscCall(PetscFree(indices)); 14759566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isnumb)); 14769566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ispart)); 147792bb6962SLisandro Dalcin } 14783ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 147992bb6962SLisandro Dalcin } 148092bb6962SLisandro Dalcin 148192bb6962SLisandro Dalcin /*@C 148292bb6962SLisandro Dalcin PCASMDestroySubdomains - Destroys the index sets created with 1483f1580f4eSBarry Smith `PCASMCreateSubdomains()`. Should be called after setting subdomains with `PCASMSetLocalSubdomains()`. 148492bb6962SLisandro Dalcin 148592bb6962SLisandro Dalcin Collective 148692bb6962SLisandro Dalcin 148792bb6962SLisandro Dalcin Input Parameters: 148892bb6962SLisandro Dalcin + n - the number of index sets 14892b691e39SMatthew Knepley . is - the array of index sets 14902fe279fdSBarry Smith - is_local - the array of local index sets, can be `NULL` 149192bb6962SLisandro Dalcin 149292bb6962SLisandro Dalcin Level: advanced 149392bb6962SLisandro Dalcin 1494562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMCreateSubdomains()`, `PCASMSetLocalSubdomains()` 149592bb6962SLisandro Dalcin @*/ 1496d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMDestroySubdomains(PetscInt n, IS is[], IS is_local[]) 1497d71ae5a4SJacob Faibussowitsch { 149892bb6962SLisandro Dalcin PetscInt i; 14995fd66863SKarl Rupp 150092bb6962SLisandro Dalcin PetscFunctionBegin; 15013ba16761SJacob Faibussowitsch if (n <= 0) PetscFunctionReturn(PETSC_SUCCESS); 1502a35b7b57SDmitry Karpeev if (is) { 15034f572ea9SToby Isaac PetscAssertPointer(is, 2); 15049566063dSJacob Faibussowitsch for (i = 0; i < n; i++) PetscCall(ISDestroy(&is[i])); 15059566063dSJacob Faibussowitsch PetscCall(PetscFree(is)); 1506a35b7b57SDmitry Karpeev } 15072b691e39SMatthew Knepley if (is_local) { 15084f572ea9SToby Isaac PetscAssertPointer(is_local, 3); 15099566063dSJacob Faibussowitsch for (i = 0; i < n; i++) PetscCall(ISDestroy(&is_local[i])); 15109566063dSJacob Faibussowitsch PetscCall(PetscFree(is_local)); 15112b691e39SMatthew Knepley } 15123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 151392bb6962SLisandro Dalcin } 151492bb6962SLisandro Dalcin 15156141accfSBarry Smith /*@C 15164b9ad928SBarry Smith PCASMCreateSubdomains2D - Creates the index sets for the overlapping Schwarz 1517f1580f4eSBarry Smith preconditioner, `PCASM`, for a two-dimensional problem on a regular grid. 15184b9ad928SBarry Smith 15194b9ad928SBarry Smith Not Collective 15204b9ad928SBarry Smith 15214b9ad928SBarry Smith Input Parameters: 15226b867d5aSJose E. Roman + m - the number of mesh points in the x direction 15236b867d5aSJose E. Roman . n - the number of mesh points in the y direction 15246b867d5aSJose E. Roman . M - the number of subdomains in the x direction 15256b867d5aSJose E. Roman . N - the number of subdomains in the y direction 15264b9ad928SBarry Smith . dof - degrees of freedom per node 15274b9ad928SBarry Smith - overlap - overlap in mesh lines 15284b9ad928SBarry Smith 15294b9ad928SBarry Smith Output Parameters: 15304b9ad928SBarry Smith + Nsub - the number of subdomains created 15313d03bb48SJed Brown . is - array of index sets defining overlapping (if overlap > 0) subdomains 15323d03bb48SJed Brown - is_local - array of index sets defining non-overlapping subdomains 15334b9ad928SBarry Smith 153420f4b53cSBarry Smith Level: advanced 153520f4b53cSBarry Smith 15364b9ad928SBarry Smith Note: 1537f1580f4eSBarry Smith Presently `PCAMSCreateSubdomains2d()` is valid only for sequential 15384b9ad928SBarry Smith preconditioners. More general related routines are 1539f1580f4eSBarry Smith `PCASMSetTotalSubdomains()` and `PCASMSetLocalSubdomains()`. 15404b9ad928SBarry Smith 1541feefa0e1SJacob Faibussowitsch Fortran Notes: 15426141accfSBarry Smith The `IS` must be declared as an array of length long enough to hold `Nsub` entries 15436141accfSBarry Smith 1544562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMSetTotalSubdomains()`, `PCASMSetLocalSubdomains()`, `PCASMGetSubKSP()`, 1545db781477SPatrick Sanan `PCASMSetOverlap()` 15464b9ad928SBarry Smith @*/ 1547d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMCreateSubdomains2D(PetscInt m, PetscInt n, PetscInt M, PetscInt N, PetscInt dof, PetscInt overlap, PetscInt *Nsub, IS **is, IS **is_local) 1548d71ae5a4SJacob Faibussowitsch { 15493d03bb48SJed Brown PetscInt i, j, height, width, ystart, xstart, yleft, yright, xleft, xright, loc_outer; 155013f74950SBarry Smith PetscInt nidx, *idx, loc, ii, jj, count; 15514b9ad928SBarry Smith 15524b9ad928SBarry Smith PetscFunctionBegin; 15537827d75bSBarry Smith PetscCheck(dof == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "dof must be 1"); 15544b9ad928SBarry Smith 15554b9ad928SBarry Smith *Nsub = N * M; 15569566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(*Nsub, is)); 15579566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(*Nsub, is_local)); 15584b9ad928SBarry Smith ystart = 0; 15593d03bb48SJed Brown loc_outer = 0; 15604b9ad928SBarry Smith for (i = 0; i < N; i++) { 15614b9ad928SBarry Smith height = n / N + ((n % N) > i); /* height of subdomain */ 15627827d75bSBarry Smith PetscCheck(height >= 2, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many N subdomains for mesh dimension n"); 15639371c9d4SSatish Balay yleft = ystart - overlap; 15649371c9d4SSatish Balay if (yleft < 0) yleft = 0; 15659371c9d4SSatish Balay yright = ystart + height + overlap; 15669371c9d4SSatish Balay if (yright > n) yright = n; 15674b9ad928SBarry Smith xstart = 0; 15684b9ad928SBarry Smith for (j = 0; j < M; j++) { 15694b9ad928SBarry Smith width = m / M + ((m % M) > j); /* width of subdomain */ 15707827d75bSBarry Smith PetscCheck(width >= 2, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many M subdomains for mesh dimension m"); 15719371c9d4SSatish Balay xleft = xstart - overlap; 15729371c9d4SSatish Balay if (xleft < 0) xleft = 0; 15739371c9d4SSatish Balay xright = xstart + width + overlap; 15749371c9d4SSatish Balay if (xright > m) xright = m; 15754b9ad928SBarry Smith nidx = (xright - xleft) * (yright - yleft); 15769566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nidx, &idx)); 15774b9ad928SBarry Smith loc = 0; 15784b9ad928SBarry Smith for (ii = yleft; ii < yright; ii++) { 15794b9ad928SBarry Smith count = m * ii + xleft; 15802fa5cd67SKarl Rupp for (jj = xleft; jj < xright; jj++) idx[loc++] = count++; 15814b9ad928SBarry Smith } 15829566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, nidx, idx, PETSC_COPY_VALUES, (*is) + loc_outer)); 15833d03bb48SJed Brown if (overlap == 0) { 15849566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)(*is)[loc_outer])); 15852fa5cd67SKarl Rupp 15863d03bb48SJed Brown (*is_local)[loc_outer] = (*is)[loc_outer]; 15873d03bb48SJed Brown } else { 15883d03bb48SJed Brown for (loc = 0, ii = ystart; ii < ystart + height; ii++) { 1589ad540459SPierre Jolivet for (jj = xstart; jj < xstart + width; jj++) idx[loc++] = m * ii + jj; 15903d03bb48SJed Brown } 15919566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, loc, idx, PETSC_COPY_VALUES, *is_local + loc_outer)); 15923d03bb48SJed Brown } 15939566063dSJacob Faibussowitsch PetscCall(PetscFree(idx)); 15944b9ad928SBarry Smith xstart += width; 15953d03bb48SJed Brown loc_outer++; 15964b9ad928SBarry Smith } 15974b9ad928SBarry Smith ystart += height; 15984b9ad928SBarry Smith } 15999566063dSJacob Faibussowitsch for (i = 0; i < *Nsub; i++) PetscCall(ISSort((*is)[i])); 16003ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 16014b9ad928SBarry Smith } 16024b9ad928SBarry Smith 16034b9ad928SBarry Smith /*@C 16044b9ad928SBarry Smith PCASMGetLocalSubdomains - Gets the local subdomains (for this processor 1605f1580f4eSBarry Smith only) for the additive Schwarz preconditioner, `PCASM`. 16064b9ad928SBarry Smith 1607ad4df100SBarry Smith Not Collective 16084b9ad928SBarry Smith 16094b9ad928SBarry Smith Input Parameter: 16104b9ad928SBarry Smith . pc - the preconditioner context 16114b9ad928SBarry Smith 16124b9ad928SBarry Smith Output Parameters: 1613f17a6784SPierre Jolivet + n - if requested, the number of subdomains for this processor (default value = 1) 1614f17a6784SPierre Jolivet . is - if requested, the index sets that define the subdomains for this processor 161520f4b53cSBarry Smith - is_local - if requested, the index sets that define the local part of the subdomains for this processor (can be `NULL`) 161620f4b53cSBarry Smith 161720f4b53cSBarry Smith Level: advanced 16184b9ad928SBarry Smith 1619f1580f4eSBarry Smith Note: 1620f1580f4eSBarry Smith The `IS` numbering is in the parallel, global numbering of the vector. 16214b9ad928SBarry Smith 1622562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMSetTotalSubdomains()`, `PCASMSetOverlap()`, `PCASMGetSubKSP()`, 1623db781477SPatrick Sanan `PCASMCreateSubdomains2D()`, `PCASMSetLocalSubdomains()`, `PCASMGetLocalSubmatrices()` 16244b9ad928SBarry Smith @*/ 1625d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMGetLocalSubdomains(PC pc, PetscInt *n, IS *is[], IS *is_local[]) 1626d71ae5a4SJacob Faibussowitsch { 16272a808120SBarry Smith PC_ASM *osm = (PC_ASM *)pc->data; 1628ace3abfcSBarry Smith PetscBool match; 16294b9ad928SBarry Smith 16304b9ad928SBarry Smith PetscFunctionBegin; 16310700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 16324f572ea9SToby Isaac if (n) PetscAssertPointer(n, 2); 16334f572ea9SToby Isaac if (is) PetscAssertPointer(is, 3); 16344f572ea9SToby Isaac if (is_local) PetscAssertPointer(is_local, 4); 16359566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCASM, &match)); 163628b400f6SJacob Faibussowitsch PetscCheck(match, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "PC is not a PCASM"); 16374b9ad928SBarry Smith if (n) *n = osm->n_local_true; 16384b9ad928SBarry Smith if (is) *is = osm->is; 16392b691e39SMatthew Knepley if (is_local) *is_local = osm->is_local; 16403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 16414b9ad928SBarry Smith } 16424b9ad928SBarry Smith 16434b9ad928SBarry Smith /*@C 16444b9ad928SBarry Smith PCASMGetLocalSubmatrices - Gets the local submatrices (for this processor 1645f1580f4eSBarry Smith only) for the additive Schwarz preconditioner, `PCASM`. 16464b9ad928SBarry Smith 1647ad4df100SBarry Smith Not Collective 16484b9ad928SBarry Smith 16494b9ad928SBarry Smith Input Parameter: 16504b9ad928SBarry Smith . pc - the preconditioner context 16514b9ad928SBarry Smith 16524b9ad928SBarry Smith Output Parameters: 1653f17a6784SPierre Jolivet + n - if requested, the number of matrices for this processor (default value = 1) 1654f17a6784SPierre Jolivet - mat - if requested, the matrices 16554b9ad928SBarry Smith 16564b9ad928SBarry Smith Level: advanced 16574b9ad928SBarry Smith 165895452b02SPatrick Sanan Notes: 1659f1580f4eSBarry Smith Call after `PCSetUp()` (or `KSPSetUp()`) but before `PCApply()` and before `PCSetUpOnBlocks()`) 1660cf739d55SBarry Smith 1661f1580f4eSBarry Smith Usually one would use `PCSetModifySubMatrices()` to change the submatrices in building the preconditioner. 1662cf739d55SBarry Smith 1663562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMSetTotalSubdomains()`, `PCASMSetOverlap()`, `PCASMGetSubKSP()`, 1664db781477SPatrick Sanan `PCASMCreateSubdomains2D()`, `PCASMSetLocalSubdomains()`, `PCASMGetLocalSubdomains()`, `PCSetModifySubMatrices()` 16654b9ad928SBarry Smith @*/ 1666d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMGetLocalSubmatrices(PC pc, PetscInt *n, Mat *mat[]) 1667d71ae5a4SJacob Faibussowitsch { 16684b9ad928SBarry Smith PC_ASM *osm; 1669ace3abfcSBarry Smith PetscBool match; 16704b9ad928SBarry Smith 16714b9ad928SBarry Smith PetscFunctionBegin; 16720700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 16734f572ea9SToby Isaac if (n) PetscAssertPointer(n, 2); 16744f572ea9SToby Isaac if (mat) PetscAssertPointer(mat, 3); 167528b400f6SJacob Faibussowitsch PetscCheck(pc->setupcalled, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must call after KSPSetUp() or PCSetUp()."); 16769566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCASM, &match)); 167792bb6962SLisandro Dalcin if (!match) { 167892bb6962SLisandro Dalcin if (n) *n = 0; 16790298fd71SBarry Smith if (mat) *mat = NULL; 168092bb6962SLisandro Dalcin } else { 16814b9ad928SBarry Smith osm = (PC_ASM *)pc->data; 16824b9ad928SBarry Smith if (n) *n = osm->n_local_true; 16834b9ad928SBarry Smith if (mat) *mat = osm->pmat; 168492bb6962SLisandro Dalcin } 16853ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 16864b9ad928SBarry Smith } 1687d709ea83SDmitry Karpeev 1688d709ea83SDmitry Karpeev /*@ 1689f1580f4eSBarry Smith PCASMSetDMSubdomains - Indicates whether to use `DMCreateDomainDecomposition()` to define the subdomains, whenever possible. 1690f1ee410cSBarry Smith 1691d709ea83SDmitry Karpeev Logically Collective 1692d709ea83SDmitry Karpeev 1693d8d19677SJose E. Roman Input Parameters: 1694d709ea83SDmitry Karpeev + pc - the preconditioner 1695f1580f4eSBarry Smith - flg - boolean indicating whether to use subdomains defined by the `DM` 1696d709ea83SDmitry Karpeev 1697d709ea83SDmitry Karpeev Options Database Key: 1698f1580f4eSBarry Smith . -pc_asm_dm_subdomains <bool> - use subdomains defined by the `DM` 1699d709ea83SDmitry Karpeev 1700d709ea83SDmitry Karpeev Level: intermediate 1701d709ea83SDmitry Karpeev 1702f1580f4eSBarry Smith Note: 1703f1580f4eSBarry Smith `PCASMSetTotalSubdomains()` and `PCASMSetOverlap()` take precedence over `PCASMSetDMSubdomains()`, 1704d709ea83SDmitry Karpeev so setting either of the first two effectively turns the latter off. 1705d709ea83SDmitry Karpeev 1706562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMGetDMSubdomains()`, `PCASMSetTotalSubdomains()`, `PCASMSetOverlap()` 1707db781477SPatrick Sanan `PCASMCreateSubdomains2D()`, `PCASMSetLocalSubdomains()`, `PCASMGetLocalSubdomains()` 1708d709ea83SDmitry Karpeev @*/ 1709d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMSetDMSubdomains(PC pc, PetscBool flg) 1710d71ae5a4SJacob Faibussowitsch { 1711d709ea83SDmitry Karpeev PC_ASM *osm = (PC_ASM *)pc->data; 1712d709ea83SDmitry Karpeev PetscBool match; 1713d709ea83SDmitry Karpeev 1714d709ea83SDmitry Karpeev PetscFunctionBegin; 1715d709ea83SDmitry Karpeev PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1716d709ea83SDmitry Karpeev PetscValidLogicalCollectiveBool(pc, flg, 2); 171728b400f6SJacob Faibussowitsch PetscCheck(!pc->setupcalled, ((PetscObject)pc)->comm, PETSC_ERR_ARG_WRONGSTATE, "Not for a setup PC."); 17189566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCASM, &match)); 1719ad540459SPierre Jolivet if (match) osm->dm_subdomains = flg; 17203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1721d709ea83SDmitry Karpeev } 1722d709ea83SDmitry Karpeev 1723d709ea83SDmitry Karpeev /*@ 1724f1580f4eSBarry Smith PCASMGetDMSubdomains - Returns flag indicating whether to use `DMCreateDomainDecomposition()` to define the subdomains, whenever possible. 1725f1580f4eSBarry Smith 1726d709ea83SDmitry Karpeev Not Collective 1727d709ea83SDmitry Karpeev 1728d709ea83SDmitry Karpeev Input Parameter: 1729d709ea83SDmitry Karpeev . pc - the preconditioner 1730d709ea83SDmitry Karpeev 1731d709ea83SDmitry Karpeev Output Parameter: 173220f4b53cSBarry Smith . flg - boolean indicating whether to use subdomains defined by the `DM` 1733d709ea83SDmitry Karpeev 1734d709ea83SDmitry Karpeev Level: intermediate 1735d709ea83SDmitry Karpeev 1736562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMSetDMSubdomains()`, `PCASMSetTotalSubdomains()`, `PCASMSetOverlap()` 1737db781477SPatrick Sanan `PCASMCreateSubdomains2D()`, `PCASMSetLocalSubdomains()`, `PCASMGetLocalSubdomains()` 1738d709ea83SDmitry Karpeev @*/ 1739d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMGetDMSubdomains(PC pc, PetscBool *flg) 1740d71ae5a4SJacob Faibussowitsch { 1741d709ea83SDmitry Karpeev PC_ASM *osm = (PC_ASM *)pc->data; 1742d709ea83SDmitry Karpeev PetscBool match; 1743d709ea83SDmitry Karpeev 1744d709ea83SDmitry Karpeev PetscFunctionBegin; 1745d709ea83SDmitry Karpeev PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 17464f572ea9SToby Isaac PetscAssertPointer(flg, 2); 17479566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCASM, &match)); 174856ea09ceSStefano Zampini if (match) *flg = osm->dm_subdomains; 174956ea09ceSStefano Zampini else *flg = PETSC_FALSE; 17503ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1751d709ea83SDmitry Karpeev } 175280ec0b7dSPatrick Sanan 1753feefa0e1SJacob Faibussowitsch /*@C 1754f1580f4eSBarry Smith PCASMGetSubMatType - Gets the matrix type used for `PCASM` subsolves, as a string. 175580ec0b7dSPatrick Sanan 175680ec0b7dSPatrick Sanan Not Collective 175780ec0b7dSPatrick Sanan 175880ec0b7dSPatrick Sanan Input Parameter: 1759f1580f4eSBarry Smith . pc - the `PC` 176080ec0b7dSPatrick Sanan 176180ec0b7dSPatrick Sanan Output Parameter: 1762feefa0e1SJacob Faibussowitsch . sub_mat_type - name of matrix type 176380ec0b7dSPatrick Sanan 176480ec0b7dSPatrick Sanan Level: advanced 176580ec0b7dSPatrick Sanan 1766562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMSetSubMatType()`, `PCSetType()`, `VecSetType()`, `MatType`, `Mat` 176780ec0b7dSPatrick Sanan @*/ 1768d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMGetSubMatType(PC pc, MatType *sub_mat_type) 1769d71ae5a4SJacob Faibussowitsch { 177056ea09ceSStefano Zampini PetscFunctionBegin; 177156ea09ceSStefano Zampini PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1772cac4c232SBarry Smith PetscTryMethod(pc, "PCASMGetSubMatType_C", (PC, MatType *), (pc, sub_mat_type)); 17733ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 177480ec0b7dSPatrick Sanan } 177580ec0b7dSPatrick Sanan 1776feefa0e1SJacob Faibussowitsch /*@C 1777f1580f4eSBarry Smith PCASMSetSubMatType - Set the type of matrix used for `PCASM` subsolves 177880ec0b7dSPatrick Sanan 1779c3339decSBarry Smith Collective 178080ec0b7dSPatrick Sanan 178180ec0b7dSPatrick Sanan Input Parameters: 1782f1580f4eSBarry Smith + pc - the `PC` object 1783f1580f4eSBarry Smith - sub_mat_type - the `MatType` 178480ec0b7dSPatrick Sanan 178580ec0b7dSPatrick Sanan Options Database Key: 1786f1580f4eSBarry Smith . -pc_asm_sub_mat_type <sub_mat_type> - Sets the matrix type used for subsolves, for example, seqaijviennacl. 1787f1580f4eSBarry Smith If you specify a base name like aijviennacl, the corresponding sequential type is assumed. 178880ec0b7dSPatrick Sanan 1789f1580f4eSBarry Smith Note: 17902fe279fdSBarry Smith See `MatType` for available types 179180ec0b7dSPatrick Sanan 179280ec0b7dSPatrick Sanan Level: advanced 179380ec0b7dSPatrick Sanan 1794562efe2eSBarry Smith .seealso: [](ch_ksp), `PCASM`, `PCASMGetSubMatType()`, `PCSetType()`, `VecSetType()`, `MatType`, `Mat` 179580ec0b7dSPatrick Sanan @*/ 1796d71ae5a4SJacob Faibussowitsch PetscErrorCode PCASMSetSubMatType(PC pc, MatType sub_mat_type) 1797d71ae5a4SJacob Faibussowitsch { 179856ea09ceSStefano Zampini PetscFunctionBegin; 179956ea09ceSStefano Zampini PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1800cac4c232SBarry Smith PetscTryMethod(pc, "PCASMSetSubMatType_C", (PC, MatType), (pc, sub_mat_type)); 18013ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 180280ec0b7dSPatrick Sanan } 1803