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 169371c9d4SSatish Balay static PetscErrorCode PCView_ASM(PC pc, PetscViewer viewer) { 1792bb6962SLisandro Dalcin PC_ASM *osm = (PC_ASM *)pc->data; 1813f74950SBarry Smith PetscMPIInt rank; 194d219a6aSLisandro Dalcin PetscInt i, bsz; 20ace3abfcSBarry Smith PetscBool iascii, isstring; 214b9ad928SBarry Smith PetscViewer sviewer; 22ed155784SPierre Jolivet PetscViewerFormat format; 23ed155784SPierre Jolivet const char *prefix; 244b9ad928SBarry Smith 254b9ad928SBarry Smith PetscFunctionBegin; 269566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 279566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERSTRING, &isstring)); 2832077d6dSBarry Smith if (iascii) { 293d03bb48SJed Brown char overlaps[256] = "user-defined overlap", blocks[256] = "total subdomain blocks not yet set"; 3063a3b9bcSJacob Faibussowitsch if (osm->overlap >= 0) PetscCall(PetscSNPrintf(overlaps, sizeof(overlaps), "amount of overlap = %" PetscInt_FMT, osm->overlap)); 3163a3b9bcSJacob Faibussowitsch if (osm->n > 0) PetscCall(PetscSNPrintf(blocks, sizeof(blocks), "total subdomain blocks = %" PetscInt_FMT, osm->n)); 329566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " %s, %s\n", blocks, overlaps)); 339566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " restriction/interpolation type - %s\n", PCASMTypes[osm->type])); 349566063dSJacob Faibussowitsch if (osm->dm_subdomains) PetscCall(PetscViewerASCIIPrintf(viewer, " Additive Schwarz: using DM to define subdomains\n")); 359566063dSJacob Faibussowitsch if (osm->loctype != PC_COMPOSITE_ADDITIVE) PetscCall(PetscViewerASCIIPrintf(viewer, " Additive Schwarz: local solve composition type - %s\n", PCCompositeTypes[osm->loctype])); 369566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank)); 379566063dSJacob Faibussowitsch PetscCall(PetscViewerGetFormat(viewer, &format)); 38ed155784SPierre Jolivet if (format != PETSC_VIEWER_ASCII_INFO_DETAIL) { 394d219a6aSLisandro Dalcin if (osm->ksp) { 409566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Local solver information for first block is in the following KSP and PC objects on rank 0:\n")); 419566063dSJacob Faibussowitsch PetscCall(PCGetOptionsPrefix(pc, &prefix)); 429566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Use -%sksp_view ::ascii_info_detail to display information for all blocks\n", prefix ? prefix : "")); 439566063dSJacob Faibussowitsch PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer)); 44dd400576SPatrick Sanan if (rank == 0) { 459566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 469566063dSJacob Faibussowitsch PetscCall(KSPView(osm->ksp[0], sviewer)); 479566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 484b9ad928SBarry Smith } 499566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer)); 504d219a6aSLisandro Dalcin } 514b9ad928SBarry Smith } else { 529566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 5363a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, " [%d] number of local blocks = %" PetscInt_FMT "\n", (int)rank, osm->n_local_true)); 549566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 559566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Local solver information for each block is in the following KSP and PC objects:\n")); 569566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 579566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "- - - - - - - - - - - - - - - - - -\n")); 589566063dSJacob Faibussowitsch PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer)); 59dd2fa690SBarry Smith for (i = 0; i < osm->n_local_true; i++) { 609566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(osm->is[i], &bsz)); 6163a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(sviewer, "[%d] local block number %" PetscInt_FMT ", size = %" PetscInt_FMT "\n", (int)rank, i, bsz)); 629566063dSJacob Faibussowitsch PetscCall(KSPView(osm->ksp[i], sviewer)); 639566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(sviewer, "- - - - - - - - - - - - - - - - - -\n")); 644b9ad928SBarry Smith } 659566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer)); 669566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 679566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 689566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 694b9ad928SBarry Smith } 704b9ad928SBarry Smith } else if (isstring) { 7163a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerStringSPrintf(viewer, " blocks=%" PetscInt_FMT ", overlap=%" PetscInt_FMT ", type=%s", osm->n, osm->overlap, PCASMTypes[osm->type])); 729566063dSJacob Faibussowitsch PetscCall(PetscViewerGetSubViewer(viewer, PETSC_COMM_SELF, &sviewer)); 739566063dSJacob Faibussowitsch if (osm->ksp) PetscCall(KSPView(osm->ksp[0], sviewer)); 749566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(viewer, PETSC_COMM_SELF, &sviewer)); 754b9ad928SBarry Smith } 764b9ad928SBarry Smith PetscFunctionReturn(0); 774b9ad928SBarry Smith } 784b9ad928SBarry Smith 799371c9d4SSatish Balay static PetscErrorCode PCASMPrintSubdomains(PC pc) { 8092bb6962SLisandro Dalcin PC_ASM *osm = (PC_ASM *)pc->data; 8192bb6962SLisandro Dalcin const char *prefix; 8292bb6962SLisandro Dalcin char fname[PETSC_MAX_PATH_LEN + 1]; 83643535aeSDmitry Karpeev PetscViewer viewer, sviewer; 84643535aeSDmitry Karpeev char *s; 8592bb6962SLisandro Dalcin PetscInt i, j, nidx; 8692bb6962SLisandro Dalcin const PetscInt *idx; 87643535aeSDmitry Karpeev PetscMPIInt rank, size; 88846783a0SBarry Smith 8992bb6962SLisandro Dalcin PetscFunctionBegin; 909566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size)); 919566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank)); 929566063dSJacob Faibussowitsch PetscCall(PCGetOptionsPrefix(pc, &prefix)); 939566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL, prefix, "-pc_asm_print_subdomains", fname, sizeof(fname), NULL)); 949566063dSJacob Faibussowitsch if (fname[0] == 0) PetscCall(PetscStrcpy(fname, "stdout")); 959566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIOpen(PetscObjectComm((PetscObject)pc), fname, &viewer)); 96643535aeSDmitry Karpeev for (i = 0; i < osm->n_local; i++) { 97643535aeSDmitry Karpeev if (i < osm->n_local_true) { 989566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(osm->is[i], &nidx)); 999566063dSJacob Faibussowitsch PetscCall(ISGetIndices(osm->is[i], &idx)); 100643535aeSDmitry Karpeev /* Print to a string viewer; no more than 15 characters per index plus 512 char for the header.*/ 10136a9e3b9SBarry Smith #define len 16 * (nidx + 1) + 512 1029566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(len, &s)); 1039566063dSJacob Faibussowitsch PetscCall(PetscViewerStringOpen(PETSC_COMM_SELF, s, len, &sviewer)); 10436a9e3b9SBarry Smith #undef len 10563a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerStringSPrintf(sviewer, "[%d:%d] Subdomain %" PetscInt_FMT " with overlap:\n", rank, size, i)); 10648a46eb9SPierre Jolivet for (j = 0; j < nidx; j++) PetscCall(PetscViewerStringSPrintf(sviewer, "%" PetscInt_FMT " ", idx[j])); 1079566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(osm->is[i], &idx)); 1089566063dSJacob Faibussowitsch PetscCall(PetscViewerStringSPrintf(sviewer, "\n")); 1099566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&sviewer)); 1109566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 11163a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "%s", s)); 1129566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 1139566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 1149566063dSJacob Faibussowitsch PetscCall(PetscFree(s)); 1152b691e39SMatthew Knepley if (osm->is_local) { 116643535aeSDmitry Karpeev /* Print to a string viewer; no more than 15 characters per index plus 512 char for the header.*/ 11736a9e3b9SBarry Smith #define len 16 * (nidx + 1) + 512 1189566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(len, &s)); 1199566063dSJacob Faibussowitsch PetscCall(PetscViewerStringOpen(PETSC_COMM_SELF, s, len, &sviewer)); 12036a9e3b9SBarry Smith #undef len 12163a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerStringSPrintf(sviewer, "[%d:%d] Subdomain %" PetscInt_FMT " without overlap:\n", rank, size, i)); 1229566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(osm->is_local[i], &nidx)); 1239566063dSJacob Faibussowitsch PetscCall(ISGetIndices(osm->is_local[i], &idx)); 12448a46eb9SPierre Jolivet for (j = 0; j < nidx; j++) PetscCall(PetscViewerStringSPrintf(sviewer, "%" PetscInt_FMT " ", idx[j])); 1259566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(osm->is_local[i], &idx)); 1269566063dSJacob Faibussowitsch PetscCall(PetscViewerStringSPrintf(sviewer, "\n")); 1279566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&sviewer)); 1289566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 12963a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "%s", s)); 1309566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 1319566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 1329566063dSJacob Faibussowitsch PetscCall(PetscFree(s)); 133643535aeSDmitry Karpeev } 1342fa5cd67SKarl Rupp } else { 135643535aeSDmitry Karpeev /* Participate in collective viewer calls. */ 1369566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 1379566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 1389566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 139643535aeSDmitry Karpeev /* Assume either all ranks have is_local or none do. */ 140643535aeSDmitry Karpeev if (osm->is_local) { 1419566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 1429566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 1439566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 144643535aeSDmitry Karpeev } 145fdc48646SDmitry Karpeev } 14692bb6962SLisandro Dalcin } 1479566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 1489566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 14992bb6962SLisandro Dalcin PetscFunctionReturn(0); 15092bb6962SLisandro Dalcin } 15192bb6962SLisandro Dalcin 1529371c9d4SSatish Balay static PetscErrorCode PCSetUp_ASM(PC pc) { 1534b9ad928SBarry Smith PC_ASM *osm = (PC_ASM *)pc->data; 15457501b6eSBarry Smith PetscBool flg; 15587e86712SBarry Smith PetscInt i, m, m_local; 1564b9ad928SBarry Smith MatReuse scall = MAT_REUSE_MATRIX; 1574b9ad928SBarry Smith IS isl; 1584b9ad928SBarry Smith KSP ksp; 1594b9ad928SBarry Smith PC subpc; 1602dcb1b2aSMatthew Knepley const char *prefix, *pprefix; 16123ce1328SBarry Smith Vec vec; 1620298fd71SBarry Smith DM *domain_dm = NULL; 1634b9ad928SBarry Smith 1644b9ad928SBarry Smith PetscFunctionBegin; 1654b9ad928SBarry Smith if (!pc->setupcalled) { 166265a2120SBarry Smith PetscInt m; 16792bb6962SLisandro Dalcin 1682b837212SDmitry Karpeev /* Note: if subdomains have been set either via PCASMSetTotalSubdomains() or via PCASMSetLocalSubdomains(), osm->n_local_true will not be PETSC_DECIDE */ 1692b837212SDmitry Karpeev if (osm->n_local_true == PETSC_DECIDE) { 17069ca1f37SDmitry Karpeev /* no subdomains given */ 17165db9045SDmitry Karpeev /* try pc->dm first, if allowed */ 172d709ea83SDmitry Karpeev if (osm->dm_subdomains && pc->dm) { 173feb221f8SDmitry Karpeev PetscInt num_domains, d; 174feb221f8SDmitry Karpeev char **domain_names; 1758d4ac253SDmitry Karpeev IS *inner_domain_is, *outer_domain_is; 1769566063dSJacob Faibussowitsch PetscCall(DMCreateDomainDecomposition(pc->dm, &num_domains, &domain_names, &inner_domain_is, &outer_domain_is, &domain_dm)); 177704f0395SPatrick Sanan osm->overlap = -1; /* We do not want to increase the overlap of the IS. 178704f0395SPatrick Sanan A future improvement of this code might allow one to use 179704f0395SPatrick Sanan DM-defined subdomains and also increase the overlap, 180704f0395SPatrick Sanan but that is not currently supported */ 1811baa6e33SBarry Smith if (num_domains) PetscCall(PCASMSetLocalSubdomains(pc, num_domains, outer_domain_is, inner_domain_is)); 182feb221f8SDmitry Karpeev for (d = 0; d < num_domains; ++d) { 1839566063dSJacob Faibussowitsch if (domain_names) PetscCall(PetscFree(domain_names[d])); 1849566063dSJacob Faibussowitsch if (inner_domain_is) PetscCall(ISDestroy(&inner_domain_is[d])); 1859566063dSJacob Faibussowitsch if (outer_domain_is) PetscCall(ISDestroy(&outer_domain_is[d])); 186feb221f8SDmitry Karpeev } 1879566063dSJacob Faibussowitsch PetscCall(PetscFree(domain_names)); 1889566063dSJacob Faibussowitsch PetscCall(PetscFree(inner_domain_is)); 1899566063dSJacob Faibussowitsch PetscCall(PetscFree(outer_domain_is)); 190feb221f8SDmitry Karpeev } 1912b837212SDmitry Karpeev if (osm->n_local_true == PETSC_DECIDE) { 19269ca1f37SDmitry Karpeev /* still no subdomains; use one subdomain per processor */ 1932b837212SDmitry Karpeev osm->n_local_true = 1; 194feb221f8SDmitry Karpeev } 1952b837212SDmitry Karpeev } 1962b837212SDmitry Karpeev { /* determine the global and max number of subdomains */ 1979371c9d4SSatish Balay struct { 1989371c9d4SSatish Balay PetscInt max, sum; 1999371c9d4SSatish Balay } inwork, outwork; 200c10200c1SHong Zhang PetscMPIInt size; 201c10200c1SHong Zhang 2026ac3741eSJed Brown inwork.max = osm->n_local_true; 2036ac3741eSJed Brown inwork.sum = osm->n_local_true; 2041c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(&inwork, &outwork, 1, MPIU_2INT, MPIU_MAXSUM_OP, PetscObjectComm((PetscObject)pc))); 2056ac3741eSJed Brown osm->n_local = outwork.max; 2066ac3741eSJed Brown osm->n = outwork.sum; 207c10200c1SHong Zhang 2089566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size)); 209c10200c1SHong Zhang if (outwork.max == 1 && outwork.sum == size) { 2107dae84e0SHong Zhang /* osm->n_local_true = 1 on all processes, set this option may enable use of optimized MatCreateSubMatrices() implementation */ 2119566063dSJacob Faibussowitsch PetscCall(MatSetOption(pc->pmat, MAT_SUBMAT_SINGLEIS, PETSC_TRUE)); 212c10200c1SHong Zhang } 2134b9ad928SBarry Smith } 21478904715SLisandro Dalcin if (!osm->is) { /* create the index sets */ 2159566063dSJacob Faibussowitsch PetscCall(PCASMCreateSubdomains(pc->pmat, osm->n_local_true, &osm->is)); 2164b9ad928SBarry Smith } 217f5234e35SJed Brown if (osm->n_local_true > 1 && !osm->is_local) { 2189566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(osm->n_local_true, &osm->is_local)); 219f5234e35SJed Brown for (i = 0; i < osm->n_local_true; i++) { 220f5234e35SJed Brown if (osm->overlap > 0) { /* With positive overlap, osm->is[i] will be modified */ 2219566063dSJacob Faibussowitsch PetscCall(ISDuplicate(osm->is[i], &osm->is_local[i])); 2229566063dSJacob Faibussowitsch PetscCall(ISCopy(osm->is[i], osm->is_local[i])); 223f5234e35SJed Brown } else { 2249566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)osm->is[i])); 225f5234e35SJed Brown osm->is_local[i] = osm->is[i]; 226f5234e35SJed Brown } 227f5234e35SJed Brown } 228f5234e35SJed Brown } 2299566063dSJacob Faibussowitsch PetscCall(PCGetOptionsPrefix(pc, &prefix)); 2303d03bb48SJed Brown if (osm->overlap > 0) { 2314b9ad928SBarry Smith /* Extend the "overlapping" regions by a number of steps */ 2329566063dSJacob Faibussowitsch PetscCall(MatIncreaseOverlap(pc->pmat, osm->n_local_true, osm->is, osm->overlap)); 2333d03bb48SJed Brown } 2346ed231c7SMatthew Knepley if (osm->sort_indices) { 23592bb6962SLisandro Dalcin for (i = 0; i < osm->n_local_true; i++) { 2369566063dSJacob Faibussowitsch PetscCall(ISSort(osm->is[i])); 23748a46eb9SPierre Jolivet if (osm->is_local) PetscCall(ISSort(osm->is_local[i])); 2384b9ad928SBarry Smith } 2396ed231c7SMatthew Knepley } 24098d3e228SPierre Jolivet flg = PETSC_FALSE; 24198d3e228SPierre Jolivet PetscCall(PetscOptionsHasName(NULL, prefix, "-pc_asm_print_subdomains", &flg)); 24298d3e228SPierre Jolivet if (flg) PetscCall(PCASMPrintSubdomains(pc)); 243f6991133SBarry Smith if (!osm->ksp) { 24478904715SLisandro Dalcin /* Create the local solvers */ 2459566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(osm->n_local_true, &osm->ksp)); 24648a46eb9SPierre Jolivet if (domain_dm) PetscCall(PetscInfo(pc, "Setting up ASM subproblems using the embedded DM\n")); 24792bb6962SLisandro Dalcin for (i = 0; i < osm->n_local_true; i++) { 2489566063dSJacob Faibussowitsch PetscCall(KSPCreate(PETSC_COMM_SELF, &ksp)); 2499566063dSJacob Faibussowitsch PetscCall(KSPSetErrorIfNotConverged(ksp, pc->erroriffailure)); 2509566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)pc, (PetscObject)ksp)); 2519566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)ksp, (PetscObject)pc, 1)); 2529566063dSJacob Faibussowitsch PetscCall(KSPSetType(ksp, KSPPREONLY)); 2539566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp, &subpc)); 2549566063dSJacob Faibussowitsch PetscCall(PCGetOptionsPrefix(pc, &prefix)); 2559566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(ksp, prefix)); 2569566063dSJacob Faibussowitsch PetscCall(KSPAppendOptionsPrefix(ksp, "sub_")); 257feb221f8SDmitry Karpeev if (domain_dm) { 2589566063dSJacob Faibussowitsch PetscCall(KSPSetDM(ksp, domain_dm[i])); 2599566063dSJacob Faibussowitsch PetscCall(KSPSetDMActive(ksp, PETSC_FALSE)); 2609566063dSJacob Faibussowitsch PetscCall(DMDestroy(&domain_dm[i])); 261feb221f8SDmitry Karpeev } 2624b9ad928SBarry Smith osm->ksp[i] = ksp; 2634b9ad928SBarry Smith } 2641baa6e33SBarry Smith if (domain_dm) PetscCall(PetscFree(domain_dm)); 265f6991133SBarry Smith } 2661dd8081eSeaulisa 2679566063dSJacob Faibussowitsch PetscCall(ISConcatenate(PETSC_COMM_SELF, osm->n_local_true, osm->is, &osm->lis)); 2689566063dSJacob Faibussowitsch PetscCall(ISSortRemoveDups(osm->lis)); 2699566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(osm->lis, &m)); 2701dd8081eSeaulisa 2714b9ad928SBarry Smith scall = MAT_INITIAL_MATRIX; 2724b9ad928SBarry Smith } else { 2734b9ad928SBarry Smith /* 2744b9ad928SBarry Smith Destroy the blocks from the previous iteration 2754b9ad928SBarry Smith */ 276e09e08ccSBarry Smith if (pc->flag == DIFFERENT_NONZERO_PATTERN) { 2779566063dSJacob Faibussowitsch PetscCall(MatDestroyMatrices(osm->n_local_true, &osm->pmat)); 2784b9ad928SBarry Smith scall = MAT_INITIAL_MATRIX; 2794b9ad928SBarry Smith } 2804b9ad928SBarry Smith } 2814b9ad928SBarry Smith 282b58cb649SBarry Smith /* Destroy previous submatrices of a different type than pc->pmat since MAT_REUSE_MATRIX won't work in that case */ 28382b5ce2aSStefano Zampini if (scall == MAT_REUSE_MATRIX && osm->sub_mat_type) { 28448a46eb9SPierre Jolivet if (osm->n_local_true > 0) PetscCall(MatDestroySubMatrices(osm->n_local_true, &osm->pmat)); 285b58cb649SBarry Smith scall = MAT_INITIAL_MATRIX; 286b58cb649SBarry Smith } 287b58cb649SBarry Smith 28892bb6962SLisandro Dalcin /* 28992bb6962SLisandro Dalcin Extract out the submatrices 29092bb6962SLisandro Dalcin */ 2919566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(pc->pmat, osm->n_local_true, osm->is, osm->is, scall, &osm->pmat)); 29292bb6962SLisandro Dalcin if (scall == MAT_INITIAL_MATRIX) { 2939566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptionsPrefix((PetscObject)pc->pmat, &pprefix)); 29492bb6962SLisandro Dalcin for (i = 0; i < osm->n_local_true; i++) { 2959566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)pc, (PetscObject)osm->pmat[i])); 2969566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)osm->pmat[i], pprefix)); 29792bb6962SLisandro Dalcin } 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 3117827d75bSBarry Smith PetscCheck(!osm->is_local || (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()"); 31248a46eb9SPierre Jolivet if (osm->is_local && osm->type == PC_ASM_RESTRICT && 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 } 4074b9ad928SBarry Smith PetscFunctionReturn(0); 4084b9ad928SBarry Smith } 4094b9ad928SBarry Smith 4109371c9d4SSatish Balay static PetscErrorCode PCSetUpOnBlocks_ASM(PC pc) { 4114b9ad928SBarry Smith PC_ASM *osm = (PC_ASM *)pc->data; 41213f74950SBarry Smith PetscInt i; 413539c167fSBarry Smith KSPConvergedReason reason; 4144b9ad928SBarry Smith 4154b9ad928SBarry Smith PetscFunctionBegin; 4164b9ad928SBarry Smith for (i = 0; i < osm->n_local_true; i++) { 4179566063dSJacob Faibussowitsch PetscCall(KSPSetUp(osm->ksp[i])); 4189566063dSJacob Faibussowitsch PetscCall(KSPGetConvergedReason(osm->ksp[i], &reason)); 419ad540459SPierre Jolivet if (reason == KSP_DIVERGED_PC_FAILED) pc->failedreason = PC_SUBPC_ERROR; 4204b9ad928SBarry Smith } 4214b9ad928SBarry Smith PetscFunctionReturn(0); 4224b9ad928SBarry Smith } 4234b9ad928SBarry Smith 4249371c9d4SSatish Balay static PetscErrorCode PCApply_ASM(PC pc, Vec x, Vec y) { 4254b9ad928SBarry Smith PC_ASM *osm = (PC_ASM *)pc->data; 4261dd8081eSeaulisa PetscInt i, n_local_true = osm->n_local_true; 4274b9ad928SBarry Smith ScatterMode forward = SCATTER_FORWARD, reverse = SCATTER_REVERSE; 4284b9ad928SBarry Smith 4294b9ad928SBarry Smith PetscFunctionBegin; 4304b9ad928SBarry Smith /* 43148e38a8aSPierre Jolivet support for limiting the restriction or interpolation to only local 4324b9ad928SBarry Smith subdomain values (leaving the other values 0). 4334b9ad928SBarry Smith */ 4344b9ad928SBarry Smith if (!(osm->type & PC_ASM_RESTRICT)) { 4354b9ad928SBarry Smith forward = SCATTER_FORWARD_LOCAL; 4364b9ad928SBarry Smith /* have to zero the work RHS since scatter may leave some slots empty */ 4379566063dSJacob Faibussowitsch PetscCall(VecSet(osm->lx, 0.0)); 4384b9ad928SBarry Smith } 439ad540459SPierre Jolivet if (!(osm->type & PC_ASM_INTERPOLATE)) reverse = SCATTER_REVERSE_LOCAL; 4404b9ad928SBarry Smith 441347574c9Seaulisa if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE || osm->loctype == PC_COMPOSITE_ADDITIVE) { 442b0de9f23SBarry Smith /* zero the global and the local solutions */ 4439566063dSJacob Faibussowitsch PetscCall(VecSet(y, 0.0)); 4449566063dSJacob Faibussowitsch PetscCall(VecSet(osm->ly, 0.0)); 445347574c9Seaulisa 44648e38a8aSPierre Jolivet /* copy the global RHS to local RHS including the ghost nodes */ 4479566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->restriction, x, osm->lx, INSERT_VALUES, forward)); 4489566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->restriction, x, osm->lx, INSERT_VALUES, forward)); 449347574c9Seaulisa 45048e38a8aSPierre Jolivet /* restrict local RHS to the overlapping 0-block RHS */ 4519566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward)); 4529566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward)); 453d8b3b5e3Seaulisa 45412cd4985SMatthew G. Knepley /* do the local solves */ 45512cd4985SMatthew G. Knepley for (i = 0; i < n_local_true; ++i) { 456b0de9f23SBarry Smith /* solve the overlapping i-block */ 4579566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PC_ApplyOnBlocks, osm->ksp[i], osm->x[i], osm->y[i], 0)); 4589566063dSJacob Faibussowitsch PetscCall(KSPSolve(osm->ksp[i], osm->x[i], osm->y[i])); 4599566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(osm->ksp[i], pc, osm->y[i])); 4609566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PC_ApplyOnBlocks, osm->ksp[i], osm->x[i], osm->y[i], 0)); 461d8b3b5e3Seaulisa 462910cf402Sprj- if (osm->lprolongation) { /* interpolate the non-overlapping i-block solution to the local solution (only for restrictive additive) */ 4639566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward)); 4649566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward)); 46548e38a8aSPierre Jolivet } else { /* interpolate the overlapping i-block solution to the local solution */ 4669566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse)); 4679566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse)); 468d8b3b5e3Seaulisa } 469347574c9Seaulisa 470347574c9Seaulisa if (i < n_local_true - 1) { 47148e38a8aSPierre Jolivet /* restrict local RHS to the overlapping (i+1)-block RHS */ 4729566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->lrestriction[i + 1], osm->lx, osm->x[i + 1], INSERT_VALUES, forward)); 4739566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->lrestriction[i + 1], osm->lx, osm->x[i + 1], INSERT_VALUES, forward)); 474347574c9Seaulisa 475347574c9Seaulisa if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE) { 4768966356dSPierre Jolivet /* update the overlapping (i+1)-block RHS using the current local solution */ 4779566063dSJacob Faibussowitsch PetscCall(MatMult(osm->lmats[i + 1], osm->ly, osm->y[i + 1])); 4789566063dSJacob Faibussowitsch PetscCall(VecAXPBY(osm->x[i + 1], -1., 1., osm->y[i + 1])); 4797c3d802fSMatthew G. Knepley } 48012cd4985SMatthew G. Knepley } 48112cd4985SMatthew G. Knepley } 48248e38a8aSPierre Jolivet /* add the local solution to the global solution including the ghost nodes */ 4839566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->restriction, osm->ly, y, ADD_VALUES, reverse)); 4849566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->restriction, osm->ly, y, ADD_VALUES, reverse)); 48598921bdaSJacob Faibussowitsch } else SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Invalid local composition type: %s", PCCompositeTypes[osm->loctype]); 4864b9ad928SBarry Smith PetscFunctionReturn(0); 4874b9ad928SBarry Smith } 4884b9ad928SBarry Smith 4899371c9d4SSatish Balay static PetscErrorCode PCMatApply_ASM(PC pc, Mat X, Mat Y) { 49048e38a8aSPierre Jolivet PC_ASM *osm = (PC_ASM *)pc->data; 49148e38a8aSPierre Jolivet Mat Z, W; 49248e38a8aSPierre Jolivet Vec x; 49348e38a8aSPierre Jolivet PetscInt i, m, N; 49448e38a8aSPierre Jolivet ScatterMode forward = SCATTER_FORWARD, reverse = SCATTER_REVERSE; 49548e38a8aSPierre Jolivet 49648e38a8aSPierre Jolivet PetscFunctionBegin; 4977827d75bSBarry Smith PetscCheck(osm->n_local_true <= 1, PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "Not yet implemented"); 49848e38a8aSPierre Jolivet /* 49948e38a8aSPierre Jolivet support for limiting the restriction or interpolation to only local 50048e38a8aSPierre Jolivet subdomain values (leaving the other values 0). 50148e38a8aSPierre Jolivet */ 50248e38a8aSPierre Jolivet if (!(osm->type & PC_ASM_RESTRICT)) { 50348e38a8aSPierre Jolivet forward = SCATTER_FORWARD_LOCAL; 50448e38a8aSPierre Jolivet /* have to zero the work RHS since scatter may leave some slots empty */ 5059566063dSJacob Faibussowitsch PetscCall(VecSet(osm->lx, 0.0)); 50648e38a8aSPierre Jolivet } 507ad540459SPierre Jolivet if (!(osm->type & PC_ASM_INTERPOLATE)) reverse = SCATTER_REVERSE_LOCAL; 5089566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(osm->x[0], &m)); 5099566063dSJacob Faibussowitsch PetscCall(MatGetSize(X, NULL, &N)); 5109566063dSJacob Faibussowitsch PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, m, N, NULL, &Z)); 51148e38a8aSPierre Jolivet if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE || osm->loctype == PC_COMPOSITE_ADDITIVE) { 51248e38a8aSPierre Jolivet /* zero the global and the local solutions */ 5139566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(Y)); 5149566063dSJacob Faibussowitsch PetscCall(VecSet(osm->ly, 0.0)); 51548e38a8aSPierre Jolivet 51648e38a8aSPierre Jolivet for (i = 0; i < N; ++i) { 5179566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumnVecRead(X, i, &x)); 51848e38a8aSPierre Jolivet /* copy the global RHS to local RHS including the ghost nodes */ 5199566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->restriction, x, osm->lx, INSERT_VALUES, forward)); 5209566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->restriction, x, osm->lx, INSERT_VALUES, forward)); 5219566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumnVecRead(X, i, &x)); 52248e38a8aSPierre Jolivet 5239566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumnVecWrite(Z, i, &x)); 52448e38a8aSPierre Jolivet /* restrict local RHS to the overlapping 0-block RHS */ 5259566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->lrestriction[0], osm->lx, x, INSERT_VALUES, forward)); 5269566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->lrestriction[0], osm->lx, x, INSERT_VALUES, forward)); 5279566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumnVecWrite(Z, i, &x)); 52848e38a8aSPierre Jolivet } 5299566063dSJacob Faibussowitsch PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, m, N, NULL, &W)); 53048e38a8aSPierre Jolivet /* solve the overlapping 0-block */ 5319566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PC_ApplyOnBlocks, osm->ksp[0], Z, W, 0)); 5329566063dSJacob Faibussowitsch PetscCall(KSPMatSolve(osm->ksp[0], Z, W)); 5339566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(osm->ksp[0], pc, NULL)); 5349566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PC_ApplyOnBlocks, osm->ksp[0], Z, W, 0)); 5359566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Z)); 53648e38a8aSPierre Jolivet 53748e38a8aSPierre Jolivet for (i = 0; i < N; ++i) { 5389566063dSJacob Faibussowitsch PetscCall(VecSet(osm->ly, 0.0)); 5399566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumnVecRead(W, i, &x)); 54048e38a8aSPierre Jolivet if (osm->lprolongation) { /* interpolate the non-overlapping 0-block solution to the local solution (only for restrictive additive) */ 5419566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->lprolongation[0], x, osm->ly, ADD_VALUES, forward)); 5429566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->lprolongation[0], x, osm->ly, ADD_VALUES, forward)); 54348e38a8aSPierre Jolivet } else { /* interpolate the overlapping 0-block solution to the local solution */ 5449566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->lrestriction[0], x, osm->ly, ADD_VALUES, reverse)); 5459566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->lrestriction[0], x, osm->ly, ADD_VALUES, reverse)); 54648e38a8aSPierre Jolivet } 5479566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumnVecRead(W, i, &x)); 54848e38a8aSPierre Jolivet 5499566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumnVecWrite(Y, i, &x)); 55048e38a8aSPierre Jolivet /* add the local solution to the global solution including the ghost nodes */ 5519566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->restriction, osm->ly, x, ADD_VALUES, reverse)); 5529566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->restriction, osm->ly, x, ADD_VALUES, reverse)); 5539566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumnVecWrite(Y, i, &x)); 55448e38a8aSPierre Jolivet } 5559566063dSJacob Faibussowitsch PetscCall(MatDestroy(&W)); 55698921bdaSJacob Faibussowitsch } else SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Invalid local composition type: %s", PCCompositeTypes[osm->loctype]); 55748e38a8aSPierre Jolivet PetscFunctionReturn(0); 55848e38a8aSPierre Jolivet } 55948e38a8aSPierre Jolivet 5609371c9d4SSatish Balay static PetscErrorCode PCApplyTranspose_ASM(PC pc, Vec x, Vec y) { 5614b9ad928SBarry Smith PC_ASM *osm = (PC_ASM *)pc->data; 5621dd8081eSeaulisa PetscInt i, n_local_true = osm->n_local_true; 5634b9ad928SBarry Smith ScatterMode forward = SCATTER_FORWARD, reverse = SCATTER_REVERSE; 5644b9ad928SBarry Smith 5654b9ad928SBarry Smith PetscFunctionBegin; 5664b9ad928SBarry Smith /* 5674b9ad928SBarry Smith Support for limiting the restriction or interpolation to only local 5684b9ad928SBarry Smith subdomain values (leaving the other values 0). 5694b9ad928SBarry Smith 5704b9ad928SBarry Smith Note: these are reversed from the PCApply_ASM() because we are applying the 5714b9ad928SBarry Smith transpose of the three terms 5724b9ad928SBarry Smith */ 573d8b3b5e3Seaulisa 5744b9ad928SBarry Smith if (!(osm->type & PC_ASM_INTERPOLATE)) { 5754b9ad928SBarry Smith forward = SCATTER_FORWARD_LOCAL; 5764b9ad928SBarry Smith /* have to zero the work RHS since scatter may leave some slots empty */ 5779566063dSJacob Faibussowitsch PetscCall(VecSet(osm->lx, 0.0)); 5784b9ad928SBarry Smith } 5792fa5cd67SKarl Rupp if (!(osm->type & PC_ASM_RESTRICT)) reverse = SCATTER_REVERSE_LOCAL; 5804b9ad928SBarry Smith 581b0de9f23SBarry Smith /* zero the global and the local solutions */ 5829566063dSJacob Faibussowitsch PetscCall(VecSet(y, 0.0)); 5839566063dSJacob Faibussowitsch PetscCall(VecSet(osm->ly, 0.0)); 584d8b3b5e3Seaulisa 585b0de9f23SBarry Smith /* Copy the global RHS to local RHS including the ghost nodes */ 5869566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->restriction, x, osm->lx, INSERT_VALUES, forward)); 5879566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->restriction, x, osm->lx, INSERT_VALUES, forward)); 588d8b3b5e3Seaulisa 589b0de9f23SBarry Smith /* Restrict local RHS to the overlapping 0-block RHS */ 5909566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward)); 5919566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward)); 592d8b3b5e3Seaulisa 5934b9ad928SBarry Smith /* do the local solves */ 594d8b3b5e3Seaulisa for (i = 0; i < n_local_true; ++i) { 595b0de9f23SBarry Smith /* solve the overlapping i-block */ 5969566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PC_ApplyOnBlocks, osm->ksp[i], osm->x[i], osm->y[i], 0)); 5979566063dSJacob Faibussowitsch PetscCall(KSPSolveTranspose(osm->ksp[i], osm->x[i], osm->y[i])); 5989566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(osm->ksp[i], pc, osm->y[i])); 5999566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PC_ApplyOnBlocks, osm->ksp[i], osm->x[i], osm->y[i], 0)); 600d8b3b5e3Seaulisa 601910cf402Sprj- if (osm->lprolongation) { /* interpolate the non-overlapping i-block solution to the local solution */ 6029566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward)); 6039566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward)); 604910cf402Sprj- } else { /* interpolate the overlapping i-block solution to the local solution */ 6059566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse)); 6069566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse)); 6074b9ad928SBarry Smith } 608d8b3b5e3Seaulisa 609d8b3b5e3Seaulisa if (i < n_local_true - 1) { 610b0de9f23SBarry Smith /* Restrict local RHS to the overlapping (i+1)-block RHS */ 6119566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->lrestriction[i + 1], osm->lx, osm->x[i + 1], INSERT_VALUES, forward)); 6129566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->lrestriction[i + 1], osm->lx, osm->x[i + 1], INSERT_VALUES, forward)); 6134b9ad928SBarry Smith } 6144b9ad928SBarry Smith } 615b0de9f23SBarry Smith /* Add the local solution to the global solution including the ghost nodes */ 6169566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->restriction, osm->ly, y, ADD_VALUES, reverse)); 6179566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->restriction, osm->ly, y, ADD_VALUES, reverse)); 6184b9ad928SBarry Smith PetscFunctionReturn(0); 6194b9ad928SBarry Smith } 6204b9ad928SBarry Smith 6219371c9d4SSatish Balay static PetscErrorCode PCReset_ASM(PC pc) { 6224b9ad928SBarry Smith PC_ASM *osm = (PC_ASM *)pc->data; 62313f74950SBarry Smith PetscInt i; 6244b9ad928SBarry Smith 6254b9ad928SBarry Smith PetscFunctionBegin; 62692bb6962SLisandro Dalcin if (osm->ksp) { 62748a46eb9SPierre Jolivet for (i = 0; i < osm->n_local_true; i++) PetscCall(KSPReset(osm->ksp[i])); 62892bb6962SLisandro Dalcin } 629e09e08ccSBarry Smith if (osm->pmat) { 63048a46eb9SPierre Jolivet if (osm->n_local_true > 0) PetscCall(MatDestroySubMatrices(osm->n_local_true, &osm->pmat)); 63192bb6962SLisandro Dalcin } 6321dd8081eSeaulisa if (osm->lrestriction) { 6339566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&osm->restriction)); 6341dd8081eSeaulisa for (i = 0; i < osm->n_local_true; i++) { 6359566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&osm->lrestriction[i])); 6369566063dSJacob Faibussowitsch if (osm->lprolongation) PetscCall(VecScatterDestroy(&osm->lprolongation[i])); 6379566063dSJacob Faibussowitsch PetscCall(VecDestroy(&osm->x[i])); 6389566063dSJacob Faibussowitsch PetscCall(VecDestroy(&osm->y[i])); 6394b9ad928SBarry Smith } 6409566063dSJacob Faibussowitsch PetscCall(PetscFree(osm->lrestriction)); 6419566063dSJacob Faibussowitsch if (osm->lprolongation) PetscCall(PetscFree(osm->lprolongation)); 6429566063dSJacob Faibussowitsch PetscCall(PetscFree(osm->x)); 6439566063dSJacob Faibussowitsch PetscCall(PetscFree(osm->y)); 64492bb6962SLisandro Dalcin } 6459566063dSJacob Faibussowitsch PetscCall(PCASMDestroySubdomains(osm->n_local_true, osm->is, osm->is_local)); 6469566063dSJacob Faibussowitsch PetscCall(ISDestroy(&osm->lis)); 6479566063dSJacob Faibussowitsch PetscCall(VecDestroy(&osm->lx)); 6489566063dSJacob Faibussowitsch PetscCall(VecDestroy(&osm->ly)); 64948a46eb9SPierre Jolivet if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE) PetscCall(MatDestroyMatrices(osm->n_local_true, &osm->lmats)); 6502fa5cd67SKarl Rupp 6519566063dSJacob Faibussowitsch PetscCall(PetscFree(osm->sub_mat_type)); 65280ec0b7dSPatrick Sanan 6530a545947SLisandro Dalcin osm->is = NULL; 6540a545947SLisandro Dalcin osm->is_local = NULL; 655e91c6855SBarry Smith PetscFunctionReturn(0); 656e91c6855SBarry Smith } 657e91c6855SBarry Smith 6589371c9d4SSatish Balay static PetscErrorCode PCDestroy_ASM(PC pc) { 659e91c6855SBarry Smith PC_ASM *osm = (PC_ASM *)pc->data; 660e91c6855SBarry Smith PetscInt i; 661e91c6855SBarry Smith 662e91c6855SBarry Smith PetscFunctionBegin; 6639566063dSJacob Faibussowitsch PetscCall(PCReset_ASM(pc)); 664e91c6855SBarry Smith if (osm->ksp) { 66548a46eb9SPierre Jolivet for (i = 0; i < osm->n_local_true; i++) PetscCall(KSPDestroy(&osm->ksp[i])); 6669566063dSJacob Faibussowitsch PetscCall(PetscFree(osm->ksp)); 667e91c6855SBarry Smith } 6689566063dSJacob Faibussowitsch PetscCall(PetscFree(pc->data)); 66996322394SPierre Jolivet 6709566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetLocalSubdomains_C", NULL)); 6719566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetTotalSubdomains_C", NULL)); 6729566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetOverlap_C", NULL)); 6739566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetType_C", NULL)); 6749566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMGetType_C", NULL)); 6759566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetLocalType_C", NULL)); 6769566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMGetLocalType_C", NULL)); 6779566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetSortIndices_C", NULL)); 6789566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMGetSubKSP_C", NULL)); 6799566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMGetSubMatType_C", NULL)); 6809566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetSubMatType_C", NULL)); 6814b9ad928SBarry Smith PetscFunctionReturn(0); 6824b9ad928SBarry Smith } 6834b9ad928SBarry Smith 6849371c9d4SSatish Balay static PetscErrorCode PCSetFromOptions_ASM(PC pc, PetscOptionItems *PetscOptionsObject) { 6854b9ad928SBarry Smith PC_ASM *osm = (PC_ASM *)pc->data; 6869dcbbd2bSBarry Smith PetscInt blocks, ovl; 68757501b6eSBarry Smith PetscBool flg; 68892bb6962SLisandro Dalcin PCASMType asmtype; 68912cd4985SMatthew G. Knepley PCCompositeType loctype; 69080ec0b7dSPatrick Sanan char sub_mat_type[256]; 6914b9ad928SBarry Smith 6924b9ad928SBarry Smith PetscFunctionBegin; 693d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "Additive Schwarz options"); 6949566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_asm_dm_subdomains", "Use DMCreateDomainDecomposition() to define subdomains", "PCASMSetDMSubdomains", osm->dm_subdomains, &osm->dm_subdomains, &flg)); 6959566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_asm_blocks", "Number of subdomains", "PCASMSetTotalSubdomains", osm->n, &blocks, &flg)); 69665db9045SDmitry Karpeev if (flg) { 6979566063dSJacob Faibussowitsch PetscCall(PCASMSetTotalSubdomains(pc, blocks, NULL, NULL)); 698d709ea83SDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 69965db9045SDmitry Karpeev } 7009566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_asm_local_blocks", "Number of local subdomains", "PCASMSetLocalSubdomains", osm->n_local_true, &blocks, &flg)); 701342c94f9SMatthew G. Knepley if (flg) { 7029566063dSJacob Faibussowitsch PetscCall(PCASMSetLocalSubdomains(pc, blocks, NULL, NULL)); 703342c94f9SMatthew G. Knepley osm->dm_subdomains = PETSC_FALSE; 704342c94f9SMatthew G. Knepley } 7059566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_asm_overlap", "Number of grid points overlap", "PCASMSetOverlap", osm->overlap, &ovl, &flg)); 70665db9045SDmitry Karpeev if (flg) { 7079566063dSJacob Faibussowitsch PetscCall(PCASMSetOverlap(pc, ovl)); 708d709ea83SDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 70965db9045SDmitry Karpeev } 71090d69ab7SBarry Smith flg = PETSC_FALSE; 7119566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-pc_asm_type", "Type of restriction/extension", "PCASMSetType", PCASMTypes, (PetscEnum)osm->type, (PetscEnum *)&asmtype, &flg)); 7129566063dSJacob Faibussowitsch if (flg) PetscCall(PCASMSetType(pc, asmtype)); 71312cd4985SMatthew G. Knepley flg = PETSC_FALSE; 7149566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-pc_asm_local_type", "Type of local solver composition", "PCASMSetLocalType", PCCompositeTypes, (PetscEnum)osm->loctype, (PetscEnum *)&loctype, &flg)); 7159566063dSJacob Faibussowitsch if (flg) PetscCall(PCASMSetLocalType(pc, loctype)); 7169566063dSJacob Faibussowitsch PetscCall(PetscOptionsFList("-pc_asm_sub_mat_type", "Subsolve Matrix Type", "PCASMSetSubMatType", MatList, NULL, sub_mat_type, 256, &flg)); 7171baa6e33SBarry Smith if (flg) PetscCall(PCASMSetSubMatType(pc, sub_mat_type)); 718d0609cedSBarry Smith PetscOptionsHeadEnd(); 7194b9ad928SBarry Smith PetscFunctionReturn(0); 7204b9ad928SBarry Smith } 7214b9ad928SBarry Smith 7229371c9d4SSatish Balay static PetscErrorCode PCASMSetLocalSubdomains_ASM(PC pc, PetscInt n, IS is[], IS is_local[]) { 7234b9ad928SBarry Smith PC_ASM *osm = (PC_ASM *)pc->data; 72492bb6962SLisandro Dalcin PetscInt i; 7254b9ad928SBarry Smith 7264b9ad928SBarry Smith PetscFunctionBegin; 72763a3b9bcSJacob Faibussowitsch PetscCheck(n >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Each process must have 1 or more blocks, n = %" PetscInt_FMT, n); 7287827d75bSBarry Smith PetscCheck(!pc->setupcalled || (n == osm->n_local_true && !is), PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "PCASMSetLocalSubdomains() should be called before calling PCSetUp()."); 729e7e72b3dSBarry Smith 7304b9ad928SBarry Smith if (!pc->setupcalled) { 73192bb6962SLisandro Dalcin if (is) { 7329566063dSJacob Faibussowitsch for (i = 0; i < n; i++) PetscCall(PetscObjectReference((PetscObject)is[i])); 73392bb6962SLisandro Dalcin } 734832fc9a5SMatthew Knepley if (is_local) { 7359566063dSJacob Faibussowitsch for (i = 0; i < n; i++) PetscCall(PetscObjectReference((PetscObject)is_local[i])); 736832fc9a5SMatthew Knepley } 7379566063dSJacob Faibussowitsch PetscCall(PCASMDestroySubdomains(osm->n_local_true, osm->is, osm->is_local)); 7382fa5cd67SKarl Rupp 7394b9ad928SBarry Smith osm->n_local_true = n; 7400a545947SLisandro Dalcin osm->is = NULL; 7410a545947SLisandro Dalcin osm->is_local = NULL; 74292bb6962SLisandro Dalcin if (is) { 7439566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &osm->is)); 7442fa5cd67SKarl Rupp for (i = 0; i < n; i++) osm->is[i] = is[i]; 7453d03bb48SJed Brown /* Flag indicating that the user has set overlapping subdomains so PCASM should not increase their size. */ 7463d03bb48SJed Brown osm->overlap = -1; 74792bb6962SLisandro Dalcin } 7482b691e39SMatthew Knepley if (is_local) { 7499566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &osm->is_local)); 7502fa5cd67SKarl Rupp for (i = 0; i < n; i++) osm->is_local[i] = is_local[i]; 751a35b7b57SDmitry Karpeev if (!is) { 7529566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(osm->n_local_true, &osm->is)); 753a35b7b57SDmitry Karpeev for (i = 0; i < osm->n_local_true; i++) { 754a35b7b57SDmitry Karpeev if (osm->overlap > 0) { /* With positive overlap, osm->is[i] will be modified */ 7559566063dSJacob Faibussowitsch PetscCall(ISDuplicate(osm->is_local[i], &osm->is[i])); 7569566063dSJacob Faibussowitsch PetscCall(ISCopy(osm->is_local[i], osm->is[i])); 757a35b7b57SDmitry Karpeev } else { 7589566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)osm->is_local[i])); 759a35b7b57SDmitry Karpeev osm->is[i] = osm->is_local[i]; 760a35b7b57SDmitry Karpeev } 761a35b7b57SDmitry Karpeev } 762a35b7b57SDmitry Karpeev } 7632b691e39SMatthew Knepley } 7644b9ad928SBarry Smith } 7654b9ad928SBarry Smith PetscFunctionReturn(0); 7664b9ad928SBarry Smith } 7674b9ad928SBarry Smith 7689371c9d4SSatish Balay static PetscErrorCode PCASMSetTotalSubdomains_ASM(PC pc, PetscInt N, IS *is, IS *is_local) { 7694b9ad928SBarry Smith PC_ASM *osm = (PC_ASM *)pc->data; 77013f74950SBarry Smith PetscMPIInt rank, size; 77178904715SLisandro Dalcin PetscInt n; 7724b9ad928SBarry Smith 7734b9ad928SBarry Smith PetscFunctionBegin; 77463a3b9bcSJacob Faibussowitsch PetscCheck(N >= 1, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Number of total blocks must be > 0, N = %" PetscInt_FMT, N); 7757827d75bSBarry 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."); 7764b9ad928SBarry Smith 7774b9ad928SBarry Smith /* 778880770e9SJed Brown Split the subdomains equally among all processors 7794b9ad928SBarry Smith */ 7809566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank)); 7819566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size)); 7824b9ad928SBarry Smith n = N / size + ((N % size) > rank); 78363a3b9bcSJacob 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); 7847827d75bSBarry Smith PetscCheck(!pc->setupcalled || n == osm->n_local_true, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONGSTATE, "PCASMSetTotalSubdomains() should be called before PCSetUp()."); 7854b9ad928SBarry Smith if (!pc->setupcalled) { 7869566063dSJacob Faibussowitsch PetscCall(PCASMDestroySubdomains(osm->n_local_true, osm->is, osm->is_local)); 7872fa5cd67SKarl Rupp 7884b9ad928SBarry Smith osm->n_local_true = n; 7890a545947SLisandro Dalcin osm->is = NULL; 7900a545947SLisandro Dalcin osm->is_local = NULL; 7914b9ad928SBarry Smith } 7924b9ad928SBarry Smith PetscFunctionReturn(0); 7934b9ad928SBarry Smith } 7944b9ad928SBarry Smith 7959371c9d4SSatish Balay static PetscErrorCode PCASMSetOverlap_ASM(PC pc, PetscInt ovl) { 79692bb6962SLisandro Dalcin PC_ASM *osm = (PC_ASM *)pc->data; 7974b9ad928SBarry Smith 7984b9ad928SBarry Smith PetscFunctionBegin; 7997827d75bSBarry Smith PetscCheck(ovl >= 0, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_OUTOFRANGE, "Negative overlap value requested"); 8007827d75bSBarry Smith PetscCheck(!pc->setupcalled || ovl == osm->overlap, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "PCASMSetOverlap() should be called before PCSetUp()."); 8012fa5cd67SKarl Rupp if (!pc->setupcalled) osm->overlap = ovl; 8024b9ad928SBarry Smith PetscFunctionReturn(0); 8034b9ad928SBarry Smith } 8044b9ad928SBarry Smith 8059371c9d4SSatish Balay static PetscErrorCode PCASMSetType_ASM(PC pc, PCASMType type) { 80692bb6962SLisandro Dalcin PC_ASM *osm = (PC_ASM *)pc->data; 8074b9ad928SBarry Smith 8084b9ad928SBarry Smith PetscFunctionBegin; 8094b9ad928SBarry Smith osm->type = type; 810bf108f30SBarry Smith osm->type_set = PETSC_TRUE; 8114b9ad928SBarry Smith PetscFunctionReturn(0); 8124b9ad928SBarry Smith } 8134b9ad928SBarry Smith 8149371c9d4SSatish Balay static PetscErrorCode PCASMGetType_ASM(PC pc, PCASMType *type) { 815c60c7ad4SBarry Smith PC_ASM *osm = (PC_ASM *)pc->data; 816c60c7ad4SBarry Smith 817c60c7ad4SBarry Smith PetscFunctionBegin; 818c60c7ad4SBarry Smith *type = osm->type; 819c60c7ad4SBarry Smith PetscFunctionReturn(0); 820c60c7ad4SBarry Smith } 821c60c7ad4SBarry Smith 8229371c9d4SSatish Balay static PetscErrorCode PCASMSetLocalType_ASM(PC pc, PCCompositeType type) { 82312cd4985SMatthew G. Knepley PC_ASM *osm = (PC_ASM *)pc->data; 82412cd4985SMatthew G. Knepley 82512cd4985SMatthew G. Knepley PetscFunctionBegin; 8267827d75bSBarry 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"); 82712cd4985SMatthew G. Knepley osm->loctype = type; 82812cd4985SMatthew G. Knepley PetscFunctionReturn(0); 82912cd4985SMatthew G. Knepley } 83012cd4985SMatthew G. Knepley 8319371c9d4SSatish Balay static PetscErrorCode PCASMGetLocalType_ASM(PC pc, PCCompositeType *type) { 83212cd4985SMatthew G. Knepley PC_ASM *osm = (PC_ASM *)pc->data; 83312cd4985SMatthew G. Knepley 83412cd4985SMatthew G. Knepley PetscFunctionBegin; 83512cd4985SMatthew G. Knepley *type = osm->loctype; 83612cd4985SMatthew G. Knepley PetscFunctionReturn(0); 83712cd4985SMatthew G. Knepley } 83812cd4985SMatthew G. Knepley 8399371c9d4SSatish Balay static PetscErrorCode PCASMSetSortIndices_ASM(PC pc, PetscBool doSort) { 8406ed231c7SMatthew Knepley PC_ASM *osm = (PC_ASM *)pc->data; 8416ed231c7SMatthew Knepley 8426ed231c7SMatthew Knepley PetscFunctionBegin; 8436ed231c7SMatthew Knepley osm->sort_indices = doSort; 8446ed231c7SMatthew Knepley PetscFunctionReturn(0); 8456ed231c7SMatthew Knepley } 8466ed231c7SMatthew Knepley 8479371c9d4SSatish Balay static PetscErrorCode PCASMGetSubKSP_ASM(PC pc, PetscInt *n_local, PetscInt *first_local, KSP **ksp) { 84892bb6962SLisandro Dalcin PC_ASM *osm = (PC_ASM *)pc->data; 8494b9ad928SBarry Smith 8504b9ad928SBarry Smith PetscFunctionBegin; 8517827d75bSBarry 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"); 8524b9ad928SBarry Smith 8532fa5cd67SKarl Rupp if (n_local) *n_local = osm->n_local_true; 85492bb6962SLisandro Dalcin if (first_local) { 8559566063dSJacob Faibussowitsch PetscCallMPI(MPI_Scan(&osm->n_local_true, first_local, 1, MPIU_INT, MPI_SUM, PetscObjectComm((PetscObject)pc))); 85692bb6962SLisandro Dalcin *first_local -= osm->n_local_true; 85792bb6962SLisandro Dalcin } 858ed155784SPierre Jolivet if (ksp) *ksp = osm->ksp; 8594b9ad928SBarry Smith PetscFunctionReturn(0); 8604b9ad928SBarry Smith } 8614b9ad928SBarry Smith 8629371c9d4SSatish Balay static PetscErrorCode PCASMGetSubMatType_ASM(PC pc, MatType *sub_mat_type) { 86380ec0b7dSPatrick Sanan PC_ASM *osm = (PC_ASM *)pc->data; 86480ec0b7dSPatrick Sanan 86580ec0b7dSPatrick Sanan PetscFunctionBegin; 86680ec0b7dSPatrick Sanan PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 86780ec0b7dSPatrick Sanan PetscValidPointer(sub_mat_type, 2); 86880ec0b7dSPatrick Sanan *sub_mat_type = osm->sub_mat_type; 86980ec0b7dSPatrick Sanan PetscFunctionReturn(0); 87080ec0b7dSPatrick Sanan } 87180ec0b7dSPatrick Sanan 8729371c9d4SSatish Balay static PetscErrorCode PCASMSetSubMatType_ASM(PC pc, MatType sub_mat_type) { 87380ec0b7dSPatrick Sanan PC_ASM *osm = (PC_ASM *)pc->data; 87480ec0b7dSPatrick Sanan 87580ec0b7dSPatrick Sanan PetscFunctionBegin; 87680ec0b7dSPatrick Sanan PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 8779566063dSJacob Faibussowitsch PetscCall(PetscFree(osm->sub_mat_type)); 8789566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(sub_mat_type, (char **)&osm->sub_mat_type)); 87980ec0b7dSPatrick Sanan PetscFunctionReturn(0); 88080ec0b7dSPatrick Sanan } 88180ec0b7dSPatrick Sanan 8824b9ad928SBarry Smith /*@C 883*f1580f4eSBarry Smith PCASMSetLocalSubdomains - Sets the local subdomains (for this processor only) for the additive Schwarz preconditioner `PCASM`. 8844b9ad928SBarry Smith 885d083f849SBarry Smith Collective on pc 8864b9ad928SBarry Smith 8874b9ad928SBarry Smith Input Parameters: 8884b9ad928SBarry Smith + pc - the preconditioner context 8894b9ad928SBarry Smith . n - the number of subdomains for this processor (default value = 1) 8908c03b21aSDmitry Karpeev . is - the index set that defines the subdomains for this processor 8910298fd71SBarry Smith (or NULL for PETSc to determine subdomains) 892f1ee410cSBarry 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 893f1ee410cSBarry Smith (or NULL to not provide these) 8944b9ad928SBarry Smith 895342c94f9SMatthew G. Knepley Options Database Key: 896342c94f9SMatthew G. Knepley To set the total number of subdomain blocks rather than specify the 897342c94f9SMatthew G. Knepley index sets, use the option 898342c94f9SMatthew G. Knepley . -pc_asm_local_blocks <blks> - Sets local blocks 899342c94f9SMatthew G. Knepley 9004b9ad928SBarry Smith Notes: 901*f1580f4eSBarry Smith The `IS` numbering is in the parallel, global numbering of the vector for both is and is_local 9024b9ad928SBarry Smith 903*f1580f4eSBarry Smith By default the `PCASM` preconditioner uses 1 block per processor. 9044b9ad928SBarry Smith 905*f1580f4eSBarry Smith Use `PCASMSetTotalSubdomains()` to set the subdomains for all processors. 9064b9ad928SBarry Smith 907*f1580f4eSBarry Smith If is_local is provided and `PCASMType` is `PC_ASM_RESTRICT` then the solution only over the is_local region is interpolated 908f1ee410cSBarry Smith back to form the global solution (this is the standard restricted additive Schwarz method) 909f1ee410cSBarry Smith 910*f1580f4eSBarry 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 911f1ee410cSBarry Smith no code to handle that case. 912f1ee410cSBarry Smith 9134b9ad928SBarry Smith Level: advanced 9144b9ad928SBarry Smith 915*f1580f4eSBarry Smith .seealso: `PCASM`, `PCASMSetTotalSubdomains()`, `PCASMSetOverlap()`, `PCASMGetSubKSP()`, 916*f1580f4eSBarry Smith `PCASMCreateSubdomains2D()`, `PCASMGetLocalSubdomains()`, `PCASMType`, `PCASMSetType()`, `PCGASM` 9174b9ad928SBarry Smith @*/ 9189371c9d4SSatish Balay PetscErrorCode PCASMSetLocalSubdomains(PC pc, PetscInt n, IS is[], IS is_local[]) { 9194b9ad928SBarry Smith PetscFunctionBegin; 9200700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 921cac4c232SBarry Smith PetscTryMethod(pc, "PCASMSetLocalSubdomains_C", (PC, PetscInt, IS[], IS[]), (pc, n, is, is_local)); 9224b9ad928SBarry Smith PetscFunctionReturn(0); 9234b9ad928SBarry Smith } 9244b9ad928SBarry Smith 9254b9ad928SBarry Smith /*@C 926feb221f8SDmitry Karpeev PCASMSetTotalSubdomains - Sets the subdomains for all processors for the 927*f1580f4eSBarry Smith additive Schwarz preconditioner, `PCASM`. 9284b9ad928SBarry Smith 929*f1580f4eSBarry Smith Collective on pc, all MPI ranks must pass in the same array of `IS` 9304b9ad928SBarry Smith 9314b9ad928SBarry Smith Input Parameters: 9324b9ad928SBarry Smith + pc - the preconditioner context 933feb221f8SDmitry Karpeev . N - the number of subdomains for all processors 934feb221f8SDmitry Karpeev . is - the index sets that define the subdomains for all processors 935dfaa0c5fSBarry Smith (or NULL to ask PETSc to determine the subdomains) 9362b691e39SMatthew Knepley - is_local - the index sets that define the local part of the subdomains for this processor 937f1ee410cSBarry Smith (or NULL to not provide this information) 9384b9ad928SBarry Smith 9394b9ad928SBarry Smith Options Database Key: 9404b9ad928SBarry Smith To set the total number of subdomain blocks rather than specify the 9414b9ad928SBarry Smith index sets, use the option 9424b9ad928SBarry Smith . -pc_asm_blocks <blks> - Sets total blocks 9434b9ad928SBarry Smith 9444b9ad928SBarry Smith Notes: 945f1ee410cSBarry Smith Currently you cannot use this to set the actual subdomains with the argument is or is_local. 9464b9ad928SBarry Smith 947*f1580f4eSBarry Smith By default the `PCASM` preconditioner uses 1 block per processor. 9484b9ad928SBarry Smith 9494b9ad928SBarry Smith These index sets cannot be destroyed until after completion of the 950*f1580f4eSBarry Smith linear solves for which the `PCASM` preconditioner is being used. 9514b9ad928SBarry Smith 952*f1580f4eSBarry Smith Use `PCASMSetLocalSubdomains()` to set local subdomains. 9534b9ad928SBarry Smith 954*f1580f4eSBarry Smith The `IS` numbering is in the parallel, global numbering of the vector for both is and is_local 9551093a601SBarry Smith 9564b9ad928SBarry Smith Level: advanced 9574b9ad928SBarry Smith 958*f1580f4eSBarry Smith .seealso: `PCASM`, `PCASMSetLocalSubdomains()`, `PCASMSetOverlap()`, `PCASMGetSubKSP()`, 959*f1580f4eSBarry Smith `PCASMCreateSubdomains2D()`, `PCGASM` 9604b9ad928SBarry Smith @*/ 9619371c9d4SSatish Balay PetscErrorCode PCASMSetTotalSubdomains(PC pc, PetscInt N, IS is[], IS is_local[]) { 9624b9ad928SBarry Smith PetscFunctionBegin; 9630700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 964cac4c232SBarry Smith PetscTryMethod(pc, "PCASMSetTotalSubdomains_C", (PC, PetscInt, IS[], IS[]), (pc, N, is, is_local)); 9654b9ad928SBarry Smith PetscFunctionReturn(0); 9664b9ad928SBarry Smith } 9674b9ad928SBarry Smith 9684b9ad928SBarry Smith /*@ 9694b9ad928SBarry Smith PCASMSetOverlap - Sets the overlap between a pair of subdomains for the 970*f1580f4eSBarry Smith additive Schwarz preconditioner, `PCASM`. 9714b9ad928SBarry Smith 972d083f849SBarry Smith Logically Collective on pc 9734b9ad928SBarry Smith 9744b9ad928SBarry Smith Input Parameters: 9754b9ad928SBarry Smith + pc - the preconditioner context 9764b9ad928SBarry Smith - ovl - the amount of overlap between subdomains (ovl >= 0, default value = 1) 9774b9ad928SBarry Smith 9784b9ad928SBarry Smith Options Database Key: 9794b9ad928SBarry Smith . -pc_asm_overlap <ovl> - Sets overlap 9804b9ad928SBarry Smith 9814b9ad928SBarry Smith Notes: 982*f1580f4eSBarry Smith By default the `PCASM` preconditioner uses 1 block per processor. To use 983*f1580f4eSBarry Smith multiple blocks per perocessor, see `PCASMSetTotalSubdomains()` and 984*f1580f4eSBarry Smith `PCASMSetLocalSubdomains()` (and the option -pc_asm_blocks <blks>). 9854b9ad928SBarry Smith 9864b9ad928SBarry Smith The overlap defaults to 1, so if one desires that no additional 9874b9ad928SBarry Smith overlap be computed beyond what may have been set with a call to 988*f1580f4eSBarry Smith `PCASMSetTotalSubdomains()` or `PCASMSetLocalSubdomains()`, then ovl 9894b9ad928SBarry Smith must be set to be 0. In particular, if one does not explicitly set 9904b9ad928SBarry Smith the subdomains an application code, then all overlap would be computed 991*f1580f4eSBarry Smith internally by PETSc, and using an overlap of 0 would result in an `PCASM` 9924b9ad928SBarry Smith variant that is equivalent to the block Jacobi preconditioner. 9934b9ad928SBarry Smith 994f1ee410cSBarry Smith The default algorithm used by PETSc to increase overlap is fast, but not scalable, 995f1ee410cSBarry Smith use the option -mat_increase_overlap_scalable when the problem and number of processes is large. 996f1ee410cSBarry Smith 9974b9ad928SBarry Smith Note that one can define initial index sets with any overlap via 998*f1580f4eSBarry Smith `PCASMSetLocalSubdomains()`; the routine 999*f1580f4eSBarry Smith `PCASMSetOverlap()` merely allows PETSc to extend that overlap further 10004b9ad928SBarry Smith if desired. 10014b9ad928SBarry Smith 10024b9ad928SBarry Smith Level: intermediate 10034b9ad928SBarry Smith 1004*f1580f4eSBarry Smith .seealso: `PCASM`, `PCASMSetTotalSubdomains()`, `PCASMSetLocalSubdomains()`, `PCASMGetSubKSP()`, 1005*f1580f4eSBarry Smith `PCASMCreateSubdomains2D()`, `PCASMGetLocalSubdomains()`, `MatIncreaseOverlap()`, `PCGASM` 10064b9ad928SBarry Smith @*/ 10079371c9d4SSatish Balay PetscErrorCode PCASMSetOverlap(PC pc, PetscInt ovl) { 10084b9ad928SBarry Smith PetscFunctionBegin; 10090700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1010c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc, ovl, 2); 1011cac4c232SBarry Smith PetscTryMethod(pc, "PCASMSetOverlap_C", (PC, PetscInt), (pc, ovl)); 10124b9ad928SBarry Smith PetscFunctionReturn(0); 10134b9ad928SBarry Smith } 10144b9ad928SBarry Smith 10154b9ad928SBarry Smith /*@ 10164b9ad928SBarry Smith PCASMSetType - Sets the type of restriction and interpolation used 1017*f1580f4eSBarry Smith for local problems in the additive Schwarz method, `PCASM`. 10184b9ad928SBarry Smith 1019d083f849SBarry Smith Logically Collective on pc 10204b9ad928SBarry Smith 10214b9ad928SBarry Smith Input Parameters: 10224b9ad928SBarry Smith + pc - the preconditioner context 1023*f1580f4eSBarry Smith - type - variant of `PCASM`, one of 10244b9ad928SBarry Smith .vb 10254b9ad928SBarry Smith PC_ASM_BASIC - full interpolation and restriction 102682b5ce2aSStefano Zampini PC_ASM_RESTRICT - full restriction, local processor interpolation (default) 10274b9ad928SBarry Smith PC_ASM_INTERPOLATE - full interpolation, local processor restriction 10284b9ad928SBarry Smith PC_ASM_NONE - local processor restriction and interpolation 10294b9ad928SBarry Smith .ve 10304b9ad928SBarry Smith 10314b9ad928SBarry Smith Options Database Key: 1032*f1580f4eSBarry Smith . -pc_asm_type [basic,restrict,interpolate,none] - Sets `PCASMType` 10334b9ad928SBarry Smith 1034*f1580f4eSBarry Smith Note: 1035*f1580f4eSBarry Smith if the is_local arguments are passed to `PCASMSetLocalSubdomains()` then they are used when `PC_ASM_RESTRICT` has been selected 1036f1ee410cSBarry Smith to limit the local processor interpolation 1037f1ee410cSBarry Smith 10384b9ad928SBarry Smith Level: intermediate 10394b9ad928SBarry Smith 1040*f1580f4eSBarry Smith .seealso: `PCASM`, `PCASMSetTotalSubdomains()`, `PCASMSetTotalSubdomains()`, `PCASMGetSubKSP()`, 1041*f1580f4eSBarry Smith `PCASMCreateSubdomains2D()`, `PCASMType`, `PCASMSetLocalType()`, `PCASMGetLocalType()`, `PCGASM` 10424b9ad928SBarry Smith @*/ 10439371c9d4SSatish Balay PetscErrorCode PCASMSetType(PC pc, PCASMType type) { 10444b9ad928SBarry Smith PetscFunctionBegin; 10450700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1046c5eb9154SBarry Smith PetscValidLogicalCollectiveEnum(pc, type, 2); 1047cac4c232SBarry Smith PetscTryMethod(pc, "PCASMSetType_C", (PC, PCASMType), (pc, type)); 10484b9ad928SBarry Smith PetscFunctionReturn(0); 10494b9ad928SBarry Smith } 10504b9ad928SBarry Smith 1051c60c7ad4SBarry Smith /*@ 1052c60c7ad4SBarry Smith PCASMGetType - Gets the type of restriction and interpolation used 1053*f1580f4eSBarry Smith for local problems in the additive Schwarz method, `PCASM`. 1054c60c7ad4SBarry Smith 1055d083f849SBarry Smith Logically Collective on pc 1056c60c7ad4SBarry Smith 1057c60c7ad4SBarry Smith Input Parameter: 1058c60c7ad4SBarry Smith . pc - the preconditioner context 1059c60c7ad4SBarry Smith 1060c60c7ad4SBarry Smith Output Parameter: 1061*f1580f4eSBarry Smith . type - variant of `PCASM`, one of 1062c60c7ad4SBarry Smith 1063c60c7ad4SBarry Smith .vb 1064c60c7ad4SBarry Smith PC_ASM_BASIC - full interpolation and restriction 1065c60c7ad4SBarry Smith PC_ASM_RESTRICT - full restriction, local processor interpolation 1066c60c7ad4SBarry Smith PC_ASM_INTERPOLATE - full interpolation, local processor restriction 1067c60c7ad4SBarry Smith PC_ASM_NONE - local processor restriction and interpolation 1068c60c7ad4SBarry Smith .ve 1069c60c7ad4SBarry Smith 1070c60c7ad4SBarry Smith Options Database Key: 1071*f1580f4eSBarry Smith . -pc_asm_type [basic,restrict,interpolate,none] - Sets `PCASM` type 1072c60c7ad4SBarry Smith 1073c60c7ad4SBarry Smith Level: intermediate 1074c60c7ad4SBarry Smith 1075*f1580f4eSBarry Smith .seealso: `PCASM`, `PCASMSetTotalSubdomains()`, `PCASMSetTotalSubdomains()`, `PCASMGetSubKSP()`, `PCGASM`, 1076db781477SPatrick Sanan `PCASMCreateSubdomains2D()`, `PCASMType`, `PCASMSetType()`, `PCASMSetLocalType()`, `PCASMGetLocalType()` 1077c60c7ad4SBarry Smith @*/ 10789371c9d4SSatish Balay PetscErrorCode PCASMGetType(PC pc, PCASMType *type) { 1079c60c7ad4SBarry Smith PetscFunctionBegin; 1080c60c7ad4SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1081cac4c232SBarry Smith PetscUseMethod(pc, "PCASMGetType_C", (PC, PCASMType *), (pc, type)); 1082c60c7ad4SBarry Smith PetscFunctionReturn(0); 1083c60c7ad4SBarry Smith } 1084c60c7ad4SBarry Smith 108512cd4985SMatthew G. Knepley /*@ 1086*f1580f4eSBarry Smith PCASMSetLocalType - Sets the type of composition used for local problems in the additive Schwarz method, `PCASM`. 108712cd4985SMatthew G. Knepley 1088d083f849SBarry Smith Logically Collective on pc 108912cd4985SMatthew G. Knepley 109012cd4985SMatthew G. Knepley Input Parameters: 109112cd4985SMatthew G. Knepley + pc - the preconditioner context 109212cd4985SMatthew G. Knepley - type - type of composition, one of 109312cd4985SMatthew G. Knepley .vb 109412cd4985SMatthew G. Knepley PC_COMPOSITE_ADDITIVE - local additive combination 109512cd4985SMatthew G. Knepley PC_COMPOSITE_MULTIPLICATIVE - local multiplicative combination 109612cd4985SMatthew G. Knepley .ve 109712cd4985SMatthew G. Knepley 109812cd4985SMatthew G. Knepley Options Database Key: 109912cd4985SMatthew G. Knepley . -pc_asm_local_type [additive,multiplicative] - Sets local solver composition type 110012cd4985SMatthew G. Knepley 110112cd4985SMatthew G. Knepley Level: intermediate 110212cd4985SMatthew G. Knepley 1103*f1580f4eSBarry Smith .seealso: `PCASM`, `PCASMSetType()`, `PCASMGetType()`, `PCASMGetLocalType()`, `PCASM`, `PCASMType`, `PCASMSetType()`, `PCASMGetType()`, `PCCompositeType` 110412cd4985SMatthew G. Knepley @*/ 11059371c9d4SSatish Balay PetscErrorCode PCASMSetLocalType(PC pc, PCCompositeType type) { 110612cd4985SMatthew G. Knepley PetscFunctionBegin; 110712cd4985SMatthew G. Knepley PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 110812cd4985SMatthew G. Knepley PetscValidLogicalCollectiveEnum(pc, type, 2); 1109cac4c232SBarry Smith PetscTryMethod(pc, "PCASMSetLocalType_C", (PC, PCCompositeType), (pc, type)); 111012cd4985SMatthew G. Knepley PetscFunctionReturn(0); 111112cd4985SMatthew G. Knepley } 111212cd4985SMatthew G. Knepley 111312cd4985SMatthew G. Knepley /*@ 1114*f1580f4eSBarry Smith PCASMGetLocalType - Gets the type of composition used for local problems in the additive Schwarz method, `PCASM`. 111512cd4985SMatthew G. Knepley 1116d083f849SBarry Smith Logically Collective on pc 111712cd4985SMatthew G. Knepley 111812cd4985SMatthew G. Knepley Input Parameter: 111912cd4985SMatthew G. Knepley . pc - the preconditioner context 112012cd4985SMatthew G. Knepley 112112cd4985SMatthew G. Knepley Output Parameter: 112212cd4985SMatthew G. Knepley . type - type of composition, one of 112312cd4985SMatthew G. Knepley .vb 112412cd4985SMatthew G. Knepley PC_COMPOSITE_ADDITIVE - local additive combination 112512cd4985SMatthew G. Knepley PC_COMPOSITE_MULTIPLICATIVE - local multiplicative combination 112612cd4985SMatthew G. Knepley .ve 112712cd4985SMatthew G. Knepley 112812cd4985SMatthew G. Knepley Options Database Key: 112912cd4985SMatthew G. Knepley . -pc_asm_local_type [additive,multiplicative] - Sets local solver composition type 113012cd4985SMatthew G. Knepley 113112cd4985SMatthew G. Knepley Level: intermediate 113212cd4985SMatthew G. Knepley 1133*f1580f4eSBarry Smith .seealso: `PCASM`, `PCASMSetType()`, `PCASMGetType()`, `PCASMSetLocalType()`, `PCASMCreate()`, `PCASMType`, `PCASMSetType()`, `PCASMGetType()`, `PCCompositeType` 113412cd4985SMatthew G. Knepley @*/ 11359371c9d4SSatish Balay PetscErrorCode PCASMGetLocalType(PC pc, PCCompositeType *type) { 113612cd4985SMatthew G. Knepley PetscFunctionBegin; 113712cd4985SMatthew G. Knepley PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 113812cd4985SMatthew G. Knepley PetscValidPointer(type, 2); 1139cac4c232SBarry Smith PetscUseMethod(pc, "PCASMGetLocalType_C", (PC, PCCompositeType *), (pc, type)); 114012cd4985SMatthew G. Knepley PetscFunctionReturn(0); 114112cd4985SMatthew G. Knepley } 114212cd4985SMatthew G. Knepley 11436ed231c7SMatthew Knepley /*@ 11446ed231c7SMatthew Knepley PCASMSetSortIndices - Determines whether subdomain indices are sorted. 11456ed231c7SMatthew Knepley 1146d083f849SBarry Smith Logically Collective on pc 11476ed231c7SMatthew Knepley 11486ed231c7SMatthew Knepley Input Parameters: 11496ed231c7SMatthew Knepley + pc - the preconditioner context 11506ed231c7SMatthew Knepley - doSort - sort the subdomain indices 11516ed231c7SMatthew Knepley 11526ed231c7SMatthew Knepley Level: intermediate 11536ed231c7SMatthew Knepley 1154*f1580f4eSBarry Smith .seealso: `PCASM`, `PCASMSetLocalSubdomains()`, `PCASMSetTotalSubdomains()`, `PCASMGetSubKSP()`, 1155db781477SPatrick Sanan `PCASMCreateSubdomains2D()` 11566ed231c7SMatthew Knepley @*/ 11579371c9d4SSatish Balay PetscErrorCode PCASMSetSortIndices(PC pc, PetscBool doSort) { 11586ed231c7SMatthew Knepley PetscFunctionBegin; 11590700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1160acfcf0e5SJed Brown PetscValidLogicalCollectiveBool(pc, doSort, 2); 1161cac4c232SBarry Smith PetscTryMethod(pc, "PCASMSetSortIndices_C", (PC, PetscBool), (pc, doSort)); 11626ed231c7SMatthew Knepley PetscFunctionReturn(0); 11636ed231c7SMatthew Knepley } 11646ed231c7SMatthew Knepley 11654b9ad928SBarry Smith /*@C 1166*f1580f4eSBarry Smith PCASMGetSubKSP - Gets the local `KSP` contexts for all blocks on 11674b9ad928SBarry Smith this processor. 11684b9ad928SBarry Smith 1169d083f849SBarry Smith Collective on pc iff first_local is requested 11704b9ad928SBarry Smith 11714b9ad928SBarry Smith Input Parameter: 11724b9ad928SBarry Smith . pc - the preconditioner context 11734b9ad928SBarry Smith 11744b9ad928SBarry Smith Output Parameters: 11750298fd71SBarry Smith + n_local - the number of blocks on this processor or NULL 11760298fd71SBarry Smith . first_local - the global number of the first block on this processor or NULL, 11770298fd71SBarry Smith all processors must request or all must pass NULL 1178*f1580f4eSBarry Smith - ksp - the array of `KSP` contexts 11794b9ad928SBarry Smith 1180*f1580f4eSBarry Smith Notes: 1181*f1580f4eSBarry Smith After `PCASMGetSubKSP()` the array of `KSP`s is not to be freed. 11824b9ad928SBarry Smith 1183*f1580f4eSBarry Smith You must call `KSPSetUp()` before calling `PCASMGetSubKSP()`. 11844b9ad928SBarry Smith 1185*f1580f4eSBarry Smith Fortran Note: 1186*f1580f4eSBarry 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. 1187d29017ddSJed Brown 11884b9ad928SBarry Smith Level: advanced 11894b9ad928SBarry Smith 1190*f1580f4eSBarry Smith .seealso: `PCASM`, `PCASMSetTotalSubdomains()`, `PCASMSetTotalSubdomains()`, `PCASMSetOverlap()`, 1191db781477SPatrick Sanan `PCASMCreateSubdomains2D()`, 11924b9ad928SBarry Smith @*/ 11939371c9d4SSatish Balay PetscErrorCode PCASMGetSubKSP(PC pc, PetscInt *n_local, PetscInt *first_local, KSP *ksp[]) { 11944b9ad928SBarry Smith PetscFunctionBegin; 11950700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1196cac4c232SBarry Smith PetscUseMethod(pc, "PCASMGetSubKSP_C", (PC, PetscInt *, PetscInt *, KSP **), (pc, n_local, first_local, ksp)); 11974b9ad928SBarry Smith PetscFunctionReturn(0); 11984b9ad928SBarry Smith } 11994b9ad928SBarry Smith 12004b9ad928SBarry Smith /*MC 12013b09bd56SBarry Smith PCASM - Use the (restricted) additive Schwarz method, each block is (approximately) solved with 1202*f1580f4eSBarry Smith its own `KSP` object. 12034b9ad928SBarry Smith 12044b9ad928SBarry Smith Options Database Keys: 1205*f1580f4eSBarry Smith + -pc_asm_blocks <blks> - Sets total blocks. Defaults to one block per MPI rank. 12064b9ad928SBarry Smith . -pc_asm_overlap <ovl> - Sets overlap 1207*f1580f4eSBarry Smith . -pc_asm_type [basic,restrict,interpolate,none] - Sets `PCASMType`, default is restrict. See `PCASMSetType()` 1208*f1580f4eSBarry Smith - -pc_asm_local_type [additive, multiplicative] - Sets `PCCompositeType`, default is additive. See `PCASMSetLocalType()` 12094b9ad928SBarry Smith 12104b9ad928SBarry Smith Level: beginner 12114b9ad928SBarry Smith 1212*f1580f4eSBarry Smith Notes: 1213*f1580f4eSBarry Smith If you run with, for example, 3 blocks on 1 processor or 3 blocks on 3 processors you 1214*f1580f4eSBarry Smith will get a different convergence rate due to the default option of -pc_asm_type restrict. Use 1215*f1580f4eSBarry Smith -pc_asm_type basic to get the same convergence behavior 1216*f1580f4eSBarry Smith 1217*f1580f4eSBarry Smith Each processor can have one or more blocks, but a block cannot be shared by more 1218*f1580f4eSBarry Smith than one processor. Use `PCGASM` for subdomains shared by multiple processes. 1219*f1580f4eSBarry Smith 1220*f1580f4eSBarry Smith To set options on the solvers for each block append -sub_ to all the `KSP`, and `PC` 1221*f1580f4eSBarry Smith options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly 1222*f1580f4eSBarry Smith 1223*f1580f4eSBarry Smith To set the options on the solvers separate for each block call `PCASMGetSubKSP()` 1224*f1580f4eSBarry Smith and set the options directly on the resulting `KSP` object (you can access its `PC` with `KSPGetPC()`) 1225*f1580f4eSBarry Smith 1226c582cd25SBarry Smith References: 1227606c0280SSatish Balay + * - M Dryja, OB Widlund, An additive variant of the Schwarz alternating method for the case of many subregions 122896a0c994SBarry Smith Courant Institute, New York University Technical report 1229606c0280SSatish Balay - * - Barry Smith, Petter Bjorstad, and William Gropp, Domain Decompositions: Parallel Multilevel Methods for Elliptic Partial Differential Equations, 123096a0c994SBarry Smith Cambridge University Press. 1231c582cd25SBarry Smith 1232*f1580f4eSBarry Smith .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, `PCASMType`, `PCCompositeType`, 1233db781477SPatrick Sanan `PCBJACOBI`, `PCASMGetSubKSP()`, `PCASMSetLocalSubdomains()`, `PCASMType`, `PCASMGetType()`, `PCASMSetLocalType()`, `PCASMGetLocalType()` 1234db781477SPatrick Sanan `PCASMSetTotalSubdomains()`, `PCSetModifySubMatrices()`, `PCASMSetOverlap()`, `PCASMSetType()`, `PCCompositeType` 12354b9ad928SBarry Smith M*/ 12364b9ad928SBarry Smith 12379371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode PCCreate_ASM(PC pc) { 12384b9ad928SBarry Smith PC_ASM *osm; 12394b9ad928SBarry Smith 12404b9ad928SBarry Smith PetscFunctionBegin; 12419566063dSJacob Faibussowitsch PetscCall(PetscNewLog(pc, &osm)); 12422fa5cd67SKarl Rupp 12434b9ad928SBarry Smith osm->n = PETSC_DECIDE; 12444b9ad928SBarry Smith osm->n_local = 0; 12452b837212SDmitry Karpeev osm->n_local_true = PETSC_DECIDE; 12464b9ad928SBarry Smith osm->overlap = 1; 12470a545947SLisandro Dalcin osm->ksp = NULL; 12480a545947SLisandro Dalcin osm->restriction = NULL; 12490a545947SLisandro Dalcin osm->lprolongation = NULL; 12500a545947SLisandro Dalcin osm->lrestriction = NULL; 12510a545947SLisandro Dalcin osm->x = NULL; 12520a545947SLisandro Dalcin osm->y = NULL; 12530a545947SLisandro Dalcin osm->is = NULL; 12540a545947SLisandro Dalcin osm->is_local = NULL; 12550a545947SLisandro Dalcin osm->mat = NULL; 12560a545947SLisandro Dalcin osm->pmat = NULL; 12574b9ad928SBarry Smith osm->type = PC_ASM_RESTRICT; 125812cd4985SMatthew G. Knepley osm->loctype = PC_COMPOSITE_ADDITIVE; 12596ed231c7SMatthew Knepley osm->sort_indices = PETSC_TRUE; 1260d709ea83SDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 126180ec0b7dSPatrick Sanan osm->sub_mat_type = NULL; 12624b9ad928SBarry Smith 126392bb6962SLisandro Dalcin pc->data = (void *)osm; 12644b9ad928SBarry Smith pc->ops->apply = PCApply_ASM; 126548e38a8aSPierre Jolivet pc->ops->matapply = PCMatApply_ASM; 12664b9ad928SBarry Smith pc->ops->applytranspose = PCApplyTranspose_ASM; 12674b9ad928SBarry Smith pc->ops->setup = PCSetUp_ASM; 1268e91c6855SBarry Smith pc->ops->reset = PCReset_ASM; 12694b9ad928SBarry Smith pc->ops->destroy = PCDestroy_ASM; 12704b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_ASM; 12714b9ad928SBarry Smith pc->ops->setuponblocks = PCSetUpOnBlocks_ASM; 12724b9ad928SBarry Smith pc->ops->view = PCView_ASM; 12730a545947SLisandro Dalcin pc->ops->applyrichardson = NULL; 12744b9ad928SBarry Smith 12759566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetLocalSubdomains_C", PCASMSetLocalSubdomains_ASM)); 12769566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetTotalSubdomains_C", PCASMSetTotalSubdomains_ASM)); 12779566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetOverlap_C", PCASMSetOverlap_ASM)); 12789566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetType_C", PCASMSetType_ASM)); 12799566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMGetType_C", PCASMGetType_ASM)); 12809566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetLocalType_C", PCASMSetLocalType_ASM)); 12819566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMGetLocalType_C", PCASMGetLocalType_ASM)); 12829566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetSortIndices_C", PCASMSetSortIndices_ASM)); 12839566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMGetSubKSP_C", PCASMGetSubKSP_ASM)); 12849566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMGetSubMatType_C", PCASMGetSubMatType_ASM)); 12859566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc, "PCASMSetSubMatType_C", PCASMSetSubMatType_ASM)); 12864b9ad928SBarry Smith PetscFunctionReturn(0); 12874b9ad928SBarry Smith } 12884b9ad928SBarry Smith 128992bb6962SLisandro Dalcin /*@C 129092bb6962SLisandro Dalcin PCASMCreateSubdomains - Creates the index sets for the overlapping Schwarz 1291*f1580f4eSBarry Smith preconditioner, `PCASM`, for any problem on a general grid. 129292bb6962SLisandro Dalcin 129392bb6962SLisandro Dalcin Collective 129492bb6962SLisandro Dalcin 129592bb6962SLisandro Dalcin Input Parameters: 129692bb6962SLisandro Dalcin + A - The global matrix operator 129792bb6962SLisandro Dalcin - n - the number of local blocks 129892bb6962SLisandro Dalcin 129992bb6962SLisandro Dalcin Output Parameters: 130092bb6962SLisandro Dalcin . outis - the array of index sets defining the subdomains 130192bb6962SLisandro Dalcin 130292bb6962SLisandro Dalcin Level: advanced 130392bb6962SLisandro Dalcin 1304*f1580f4eSBarry Smith Note: 1305*f1580f4eSBarry Smith This generates nonoverlapping subdomains; the `PCASM` will generate the overlap 1306*f1580f4eSBarry Smith from these if you use `PCASMSetLocalSubdomains()` 13077d6bfa3bSBarry Smith 1308*f1580f4eSBarry Smith Fortran Note: 1309*f1580f4eSBarry Smith You must provide the array outis[] already allocated of length n. 13107d6bfa3bSBarry Smith 1311*f1580f4eSBarry Smith .seealso: `PCASM`, `PCASMSetLocalSubdomains()`, `PCASMDestroySubdomains()` 131292bb6962SLisandro Dalcin @*/ 13139371c9d4SSatish Balay PetscErrorCode PCASMCreateSubdomains(Mat A, PetscInt n, IS *outis[]) { 131492bb6962SLisandro Dalcin MatPartitioning mpart; 131592bb6962SLisandro Dalcin const char *prefix; 131692bb6962SLisandro Dalcin PetscInt i, j, rstart, rend, bs; 1317976e8c5aSLisandro Dalcin PetscBool hasop, isbaij = PETSC_FALSE, foundpart = PETSC_FALSE; 13180298fd71SBarry Smith Mat Ad = NULL, adj; 131992bb6962SLisandro Dalcin IS ispart, isnumb, *is; 132092bb6962SLisandro Dalcin 132192bb6962SLisandro Dalcin PetscFunctionBegin; 13220700a824SBarry Smith PetscValidHeaderSpecific(A, MAT_CLASSID, 1); 132392bb6962SLisandro Dalcin PetscValidPointer(outis, 3); 132463a3b9bcSJacob Faibussowitsch PetscCheck(n >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "number of local blocks must be > 0, n = %" PetscInt_FMT, n); 132592bb6962SLisandro Dalcin 132692bb6962SLisandro Dalcin /* Get prefix, row distribution, and block size */ 13279566063dSJacob Faibussowitsch PetscCall(MatGetOptionsPrefix(A, &prefix)); 13289566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A, &rstart, &rend)); 13299566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(A, &bs)); 133063a3b9bcSJacob 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); 133165e19b50SBarry Smith 133292bb6962SLisandro Dalcin /* Get diagonal block from matrix if possible */ 13339566063dSJacob Faibussowitsch PetscCall(MatHasOperation(A, MATOP_GET_DIAGONAL_BLOCK, &hasop)); 133448a46eb9SPierre Jolivet if (hasop) PetscCall(MatGetDiagonalBlock(A, &Ad)); 133592bb6962SLisandro Dalcin if (Ad) { 13369566063dSJacob Faibussowitsch PetscCall(PetscObjectBaseTypeCompare((PetscObject)Ad, MATSEQBAIJ, &isbaij)); 13379566063dSJacob Faibussowitsch if (!isbaij) PetscCall(PetscObjectBaseTypeCompare((PetscObject)Ad, MATSEQSBAIJ, &isbaij)); 133892bb6962SLisandro Dalcin } 133992bb6962SLisandro Dalcin if (Ad && n > 1) { 1340ace3abfcSBarry Smith PetscBool match, done; 134192bb6962SLisandro Dalcin /* Try to setup a good matrix partitioning if available */ 13429566063dSJacob Faibussowitsch PetscCall(MatPartitioningCreate(PETSC_COMM_SELF, &mpart)); 13439566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)mpart, prefix)); 13449566063dSJacob Faibussowitsch PetscCall(MatPartitioningSetFromOptions(mpart)); 13459566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)mpart, MATPARTITIONINGCURRENT, &match)); 134648a46eb9SPierre Jolivet if (!match) PetscCall(PetscObjectTypeCompare((PetscObject)mpart, MATPARTITIONINGSQUARE, &match)); 134792bb6962SLisandro Dalcin if (!match) { /* assume a "good" partitioner is available */ 13481a83f524SJed Brown PetscInt na; 13491a83f524SJed Brown const PetscInt *ia, *ja; 13509566063dSJacob Faibussowitsch PetscCall(MatGetRowIJ(Ad, 0, PETSC_TRUE, isbaij, &na, &ia, &ja, &done)); 135192bb6962SLisandro Dalcin if (done) { 135292bb6962SLisandro Dalcin /* Build adjacency matrix by hand. Unfortunately a call to 135392bb6962SLisandro Dalcin MatConvert(Ad,MATMPIADJ,MAT_INITIAL_MATRIX,&adj) will 135492bb6962SLisandro Dalcin remove the block-aij structure and we cannot expect 135592bb6962SLisandro Dalcin MatPartitioning to split vertices as we need */ 13560a545947SLisandro Dalcin PetscInt i, j, len, nnz, cnt, *iia = NULL, *jja = NULL; 13571a83f524SJed Brown const PetscInt *row; 135892bb6962SLisandro Dalcin nnz = 0; 135992bb6962SLisandro Dalcin for (i = 0; i < na; i++) { /* count number of nonzeros */ 136092bb6962SLisandro Dalcin len = ia[i + 1] - ia[i]; 136192bb6962SLisandro Dalcin row = ja + ia[i]; 136292bb6962SLisandro Dalcin for (j = 0; j < len; j++) { 136392bb6962SLisandro Dalcin if (row[j] == i) { /* don't count diagonal */ 13649371c9d4SSatish Balay len--; 13659371c9d4SSatish Balay break; 136692bb6962SLisandro Dalcin } 136792bb6962SLisandro Dalcin } 136892bb6962SLisandro Dalcin nnz += len; 136992bb6962SLisandro Dalcin } 13709566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(na + 1, &iia)); 13719566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nnz, &jja)); 137292bb6962SLisandro Dalcin nnz = 0; 137392bb6962SLisandro Dalcin iia[0] = 0; 137492bb6962SLisandro Dalcin for (i = 0; i < na; i++) { /* fill adjacency */ 137592bb6962SLisandro Dalcin cnt = 0; 137692bb6962SLisandro Dalcin len = ia[i + 1] - ia[i]; 137792bb6962SLisandro Dalcin row = ja + ia[i]; 137892bb6962SLisandro Dalcin for (j = 0; j < len; j++) { 137992bb6962SLisandro Dalcin if (row[j] != i) { /* if not diagonal */ 138092bb6962SLisandro Dalcin jja[nnz + cnt++] = row[j]; 138192bb6962SLisandro Dalcin } 138292bb6962SLisandro Dalcin } 138392bb6962SLisandro Dalcin nnz += cnt; 138492bb6962SLisandro Dalcin iia[i + 1] = nnz; 138592bb6962SLisandro Dalcin } 138692bb6962SLisandro Dalcin /* Partitioning of the adjacency matrix */ 13879566063dSJacob Faibussowitsch PetscCall(MatCreateMPIAdj(PETSC_COMM_SELF, na, na, iia, jja, NULL, &adj)); 13889566063dSJacob Faibussowitsch PetscCall(MatPartitioningSetAdjacency(mpart, adj)); 13899566063dSJacob Faibussowitsch PetscCall(MatPartitioningSetNParts(mpart, n)); 13909566063dSJacob Faibussowitsch PetscCall(MatPartitioningApply(mpart, &ispart)); 13919566063dSJacob Faibussowitsch PetscCall(ISPartitioningToNumbering(ispart, &isnumb)); 13929566063dSJacob Faibussowitsch PetscCall(MatDestroy(&adj)); 139392bb6962SLisandro Dalcin foundpart = PETSC_TRUE; 139492bb6962SLisandro Dalcin } 13959566063dSJacob Faibussowitsch PetscCall(MatRestoreRowIJ(Ad, 0, PETSC_TRUE, isbaij, &na, &ia, &ja, &done)); 139692bb6962SLisandro Dalcin } 13979566063dSJacob Faibussowitsch PetscCall(MatPartitioningDestroy(&mpart)); 139892bb6962SLisandro Dalcin } 139992bb6962SLisandro Dalcin 14009566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &is)); 140192bb6962SLisandro Dalcin *outis = is; 140292bb6962SLisandro Dalcin 140392bb6962SLisandro Dalcin if (!foundpart) { 140492bb6962SLisandro Dalcin /* Partitioning by contiguous chunks of rows */ 140592bb6962SLisandro Dalcin 140692bb6962SLisandro Dalcin PetscInt mbs = (rend - rstart) / bs; 140792bb6962SLisandro Dalcin PetscInt start = rstart; 140892bb6962SLisandro Dalcin for (i = 0; i < n; i++) { 140992bb6962SLisandro Dalcin PetscInt count = (mbs / n + ((mbs % n) > i)) * bs; 14109566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF, count, start, 1, &is[i])); 141192bb6962SLisandro Dalcin start += count; 141292bb6962SLisandro Dalcin } 141392bb6962SLisandro Dalcin 141492bb6962SLisandro Dalcin } else { 141592bb6962SLisandro Dalcin /* Partitioning by adjacency of diagonal block */ 141692bb6962SLisandro Dalcin 141792bb6962SLisandro Dalcin const PetscInt *numbering; 141892bb6962SLisandro Dalcin PetscInt *count, nidx, *indices, *newidx, start = 0; 141992bb6962SLisandro Dalcin /* Get node count in each partition */ 14209566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n, &count)); 14219566063dSJacob Faibussowitsch PetscCall(ISPartitioningCount(ispart, n, count)); 142292bb6962SLisandro Dalcin if (isbaij && bs > 1) { /* adjust for the block-aij case */ 142392bb6962SLisandro Dalcin for (i = 0; i < n; i++) count[i] *= bs; 142492bb6962SLisandro Dalcin } 142592bb6962SLisandro Dalcin /* Build indices from node numbering */ 14269566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(isnumb, &nidx)); 14279566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nidx, &indices)); 142892bb6962SLisandro Dalcin for (i = 0; i < nidx; i++) indices[i] = i; /* needs to be initialized */ 14299566063dSJacob Faibussowitsch PetscCall(ISGetIndices(isnumb, &numbering)); 14309566063dSJacob Faibussowitsch PetscCall(PetscSortIntWithPermutation(nidx, numbering, indices)); 14319566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(isnumb, &numbering)); 143292bb6962SLisandro Dalcin if (isbaij && bs > 1) { /* adjust for the block-aij case */ 14339566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nidx * bs, &newidx)); 14342fa5cd67SKarl Rupp for (i = 0; i < nidx; i++) { 14352fa5cd67SKarl Rupp for (j = 0; j < bs; j++) newidx[i * bs + j] = indices[i] * bs + j; 14362fa5cd67SKarl Rupp } 14379566063dSJacob Faibussowitsch PetscCall(PetscFree(indices)); 143892bb6962SLisandro Dalcin nidx *= bs; 143992bb6962SLisandro Dalcin indices = newidx; 144092bb6962SLisandro Dalcin } 144192bb6962SLisandro Dalcin /* Shift to get global indices */ 144292bb6962SLisandro Dalcin for (i = 0; i < nidx; i++) indices[i] += rstart; 144392bb6962SLisandro Dalcin 144492bb6962SLisandro Dalcin /* Build the index sets for each block */ 144592bb6962SLisandro Dalcin for (i = 0; i < n; i++) { 14469566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, count[i], &indices[start], PETSC_COPY_VALUES, &is[i])); 14479566063dSJacob Faibussowitsch PetscCall(ISSort(is[i])); 144892bb6962SLisandro Dalcin start += count[i]; 144992bb6962SLisandro Dalcin } 145092bb6962SLisandro Dalcin 14519566063dSJacob Faibussowitsch PetscCall(PetscFree(count)); 14529566063dSJacob Faibussowitsch PetscCall(PetscFree(indices)); 14539566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isnumb)); 14549566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ispart)); 145592bb6962SLisandro Dalcin } 145692bb6962SLisandro Dalcin PetscFunctionReturn(0); 145792bb6962SLisandro Dalcin } 145892bb6962SLisandro Dalcin 145992bb6962SLisandro Dalcin /*@C 146092bb6962SLisandro Dalcin PCASMDestroySubdomains - Destroys the index sets created with 1461*f1580f4eSBarry Smith `PCASMCreateSubdomains()`. Should be called after setting subdomains with `PCASMSetLocalSubdomains()`. 146292bb6962SLisandro Dalcin 146392bb6962SLisandro Dalcin Collective 146492bb6962SLisandro Dalcin 146592bb6962SLisandro Dalcin Input Parameters: 146692bb6962SLisandro Dalcin + n - the number of index sets 14672b691e39SMatthew Knepley . is - the array of index sets 14680298fd71SBarry Smith - is_local - the array of local index sets, can be NULL 146992bb6962SLisandro Dalcin 147092bb6962SLisandro Dalcin Level: advanced 147192bb6962SLisandro Dalcin 1472*f1580f4eSBarry Smith .seealso: `PCASM`, `PCASMCreateSubdomains()`, `PCASMSetLocalSubdomains()` 147392bb6962SLisandro Dalcin @*/ 14749371c9d4SSatish Balay PetscErrorCode PCASMDestroySubdomains(PetscInt n, IS is[], IS is_local[]) { 147592bb6962SLisandro Dalcin PetscInt i; 14765fd66863SKarl Rupp 147792bb6962SLisandro Dalcin PetscFunctionBegin; 1478a35b7b57SDmitry Karpeev if (n <= 0) PetscFunctionReturn(0); 1479a35b7b57SDmitry Karpeev if (is) { 148092bb6962SLisandro Dalcin PetscValidPointer(is, 2); 14819566063dSJacob Faibussowitsch for (i = 0; i < n; i++) PetscCall(ISDestroy(&is[i])); 14829566063dSJacob Faibussowitsch PetscCall(PetscFree(is)); 1483a35b7b57SDmitry Karpeev } 14842b691e39SMatthew Knepley if (is_local) { 14852b691e39SMatthew Knepley PetscValidPointer(is_local, 3); 14869566063dSJacob Faibussowitsch for (i = 0; i < n; i++) PetscCall(ISDestroy(&is_local[i])); 14879566063dSJacob Faibussowitsch PetscCall(PetscFree(is_local)); 14882b691e39SMatthew Knepley } 148992bb6962SLisandro Dalcin PetscFunctionReturn(0); 149092bb6962SLisandro Dalcin } 149192bb6962SLisandro Dalcin 14924b9ad928SBarry Smith /*@ 14934b9ad928SBarry Smith PCASMCreateSubdomains2D - Creates the index sets for the overlapping Schwarz 1494*f1580f4eSBarry Smith preconditioner, `PCASM`, for a two-dimensional problem on a regular grid. 14954b9ad928SBarry Smith 14964b9ad928SBarry Smith Not Collective 14974b9ad928SBarry Smith 14984b9ad928SBarry Smith Input Parameters: 14996b867d5aSJose E. Roman + m - the number of mesh points in the x direction 15006b867d5aSJose E. Roman . n - the number of mesh points in the y direction 15016b867d5aSJose E. Roman . M - the number of subdomains in the x direction 15026b867d5aSJose E. Roman . N - the number of subdomains in the y direction 15034b9ad928SBarry Smith . dof - degrees of freedom per node 15044b9ad928SBarry Smith - overlap - overlap in mesh lines 15054b9ad928SBarry Smith 15064b9ad928SBarry Smith Output Parameters: 15074b9ad928SBarry Smith + Nsub - the number of subdomains created 15083d03bb48SJed Brown . is - array of index sets defining overlapping (if overlap > 0) subdomains 15093d03bb48SJed Brown - is_local - array of index sets defining non-overlapping subdomains 15104b9ad928SBarry Smith 15114b9ad928SBarry Smith Note: 1512*f1580f4eSBarry Smith Presently `PCAMSCreateSubdomains2d()` is valid only for sequential 15134b9ad928SBarry Smith preconditioners. More general related routines are 1514*f1580f4eSBarry Smith `PCASMSetTotalSubdomains()` and `PCASMSetLocalSubdomains()`. 15154b9ad928SBarry Smith 15164b9ad928SBarry Smith Level: advanced 15174b9ad928SBarry Smith 1518*f1580f4eSBarry Smith .seealso: `PCASM`, `PCASMSetTotalSubdomains()`, `PCASMSetLocalSubdomains()`, `PCASMGetSubKSP()`, 1519db781477SPatrick Sanan `PCASMSetOverlap()` 15204b9ad928SBarry Smith @*/ 15219371c9d4SSatish Balay PetscErrorCode PCASMCreateSubdomains2D(PetscInt m, PetscInt n, PetscInt M, PetscInt N, PetscInt dof, PetscInt overlap, PetscInt *Nsub, IS **is, IS **is_local) { 15223d03bb48SJed Brown PetscInt i, j, height, width, ystart, xstart, yleft, yright, xleft, xright, loc_outer; 152313f74950SBarry Smith PetscInt nidx, *idx, loc, ii, jj, count; 15244b9ad928SBarry Smith 15254b9ad928SBarry Smith PetscFunctionBegin; 15267827d75bSBarry Smith PetscCheck(dof == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "dof must be 1"); 15274b9ad928SBarry Smith 15284b9ad928SBarry Smith *Nsub = N * M; 15299566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(*Nsub, is)); 15309566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(*Nsub, is_local)); 15314b9ad928SBarry Smith ystart = 0; 15323d03bb48SJed Brown loc_outer = 0; 15334b9ad928SBarry Smith for (i = 0; i < N; i++) { 15344b9ad928SBarry Smith height = n / N + ((n % N) > i); /* height of subdomain */ 15357827d75bSBarry Smith PetscCheck(height >= 2, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many N subdomains for mesh dimension n"); 15369371c9d4SSatish Balay yleft = ystart - overlap; 15379371c9d4SSatish Balay if (yleft < 0) yleft = 0; 15389371c9d4SSatish Balay yright = ystart + height + overlap; 15399371c9d4SSatish Balay if (yright > n) yright = n; 15404b9ad928SBarry Smith xstart = 0; 15414b9ad928SBarry Smith for (j = 0; j < M; j++) { 15424b9ad928SBarry Smith width = m / M + ((m % M) > j); /* width of subdomain */ 15437827d75bSBarry Smith PetscCheck(width >= 2, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Too many M subdomains for mesh dimension m"); 15449371c9d4SSatish Balay xleft = xstart - overlap; 15459371c9d4SSatish Balay if (xleft < 0) xleft = 0; 15469371c9d4SSatish Balay xright = xstart + width + overlap; 15479371c9d4SSatish Balay if (xright > m) xright = m; 15484b9ad928SBarry Smith nidx = (xright - xleft) * (yright - yleft); 15499566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nidx, &idx)); 15504b9ad928SBarry Smith loc = 0; 15514b9ad928SBarry Smith for (ii = yleft; ii < yright; ii++) { 15524b9ad928SBarry Smith count = m * ii + xleft; 15532fa5cd67SKarl Rupp for (jj = xleft; jj < xright; jj++) idx[loc++] = count++; 15544b9ad928SBarry Smith } 15559566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, nidx, idx, PETSC_COPY_VALUES, (*is) + loc_outer)); 15563d03bb48SJed Brown if (overlap == 0) { 15579566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)(*is)[loc_outer])); 15582fa5cd67SKarl Rupp 15593d03bb48SJed Brown (*is_local)[loc_outer] = (*is)[loc_outer]; 15603d03bb48SJed Brown } else { 15613d03bb48SJed Brown for (loc = 0, ii = ystart; ii < ystart + height; ii++) { 1562ad540459SPierre Jolivet for (jj = xstart; jj < xstart + width; jj++) idx[loc++] = m * ii + jj; 15633d03bb48SJed Brown } 15649566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF, loc, idx, PETSC_COPY_VALUES, *is_local + loc_outer)); 15653d03bb48SJed Brown } 15669566063dSJacob Faibussowitsch PetscCall(PetscFree(idx)); 15674b9ad928SBarry Smith xstart += width; 15683d03bb48SJed Brown loc_outer++; 15694b9ad928SBarry Smith } 15704b9ad928SBarry Smith ystart += height; 15714b9ad928SBarry Smith } 15729566063dSJacob Faibussowitsch for (i = 0; i < *Nsub; i++) PetscCall(ISSort((*is)[i])); 15734b9ad928SBarry Smith PetscFunctionReturn(0); 15744b9ad928SBarry Smith } 15754b9ad928SBarry Smith 15764b9ad928SBarry Smith /*@C 15774b9ad928SBarry Smith PCASMGetLocalSubdomains - Gets the local subdomains (for this processor 1578*f1580f4eSBarry Smith only) for the additive Schwarz preconditioner, `PCASM`. 15794b9ad928SBarry Smith 1580ad4df100SBarry Smith Not Collective 15814b9ad928SBarry Smith 15824b9ad928SBarry Smith Input Parameter: 15834b9ad928SBarry Smith . pc - the preconditioner context 15844b9ad928SBarry Smith 15854b9ad928SBarry Smith Output Parameters: 1586f17a6784SPierre Jolivet + n - if requested, the number of subdomains for this processor (default value = 1) 1587f17a6784SPierre Jolivet . is - if requested, the index sets that define the subdomains for this processor 1588f17a6784SPierre Jolivet - is_local - if requested, the index sets that define the local part of the subdomains for this processor (can be NULL) 15894b9ad928SBarry Smith 1590*f1580f4eSBarry Smith Note: 1591*f1580f4eSBarry Smith The `IS` numbering is in the parallel, global numbering of the vector. 15924b9ad928SBarry Smith 15934b9ad928SBarry Smith Level: advanced 15944b9ad928SBarry Smith 1595*f1580f4eSBarry Smith .seealso: `PCASM`, `PCASMSetTotalSubdomains()`, `PCASMSetOverlap()`, `PCASMGetSubKSP()`, 1596db781477SPatrick Sanan `PCASMCreateSubdomains2D()`, `PCASMSetLocalSubdomains()`, `PCASMGetLocalSubmatrices()` 15974b9ad928SBarry Smith @*/ 15989371c9d4SSatish Balay PetscErrorCode PCASMGetLocalSubdomains(PC pc, PetscInt *n, IS *is[], IS *is_local[]) { 15992a808120SBarry Smith PC_ASM *osm = (PC_ASM *)pc->data; 1600ace3abfcSBarry Smith PetscBool match; 16014b9ad928SBarry Smith 16024b9ad928SBarry Smith PetscFunctionBegin; 16030700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1604f17a6784SPierre Jolivet if (n) PetscValidIntPointer(n, 2); 160592bb6962SLisandro Dalcin if (is) PetscValidPointer(is, 3); 1606f17a6784SPierre Jolivet if (is_local) PetscValidPointer(is_local, 4); 16079566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCASM, &match)); 160828b400f6SJacob Faibussowitsch PetscCheck(match, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "PC is not a PCASM"); 16094b9ad928SBarry Smith if (n) *n = osm->n_local_true; 16104b9ad928SBarry Smith if (is) *is = osm->is; 16112b691e39SMatthew Knepley if (is_local) *is_local = osm->is_local; 16124b9ad928SBarry Smith PetscFunctionReturn(0); 16134b9ad928SBarry Smith } 16144b9ad928SBarry Smith 16154b9ad928SBarry Smith /*@C 16164b9ad928SBarry Smith PCASMGetLocalSubmatrices - Gets the local submatrices (for this processor 1617*f1580f4eSBarry Smith only) for the additive Schwarz preconditioner, `PCASM`. 16184b9ad928SBarry Smith 1619ad4df100SBarry Smith Not Collective 16204b9ad928SBarry Smith 16214b9ad928SBarry Smith Input Parameter: 16224b9ad928SBarry Smith . pc - the preconditioner context 16234b9ad928SBarry Smith 16244b9ad928SBarry Smith Output Parameters: 1625f17a6784SPierre Jolivet + n - if requested, the number of matrices for this processor (default value = 1) 1626f17a6784SPierre Jolivet - mat - if requested, the matrices 16274b9ad928SBarry Smith 16284b9ad928SBarry Smith Level: advanced 16294b9ad928SBarry Smith 163095452b02SPatrick Sanan Notes: 1631*f1580f4eSBarry Smith Call after `PCSetUp()` (or `KSPSetUp()`) but before `PCApply()` and before `PCSetUpOnBlocks()`) 1632cf739d55SBarry Smith 1633*f1580f4eSBarry Smith Usually one would use `PCSetModifySubMatrices()` to change the submatrices in building the preconditioner. 1634cf739d55SBarry Smith 1635*f1580f4eSBarry Smith .seealso: `PCASM`, `PCASMSetTotalSubdomains()`, `PCASMSetOverlap()`, `PCASMGetSubKSP()`, 1636db781477SPatrick Sanan `PCASMCreateSubdomains2D()`, `PCASMSetLocalSubdomains()`, `PCASMGetLocalSubdomains()`, `PCSetModifySubMatrices()` 16374b9ad928SBarry Smith @*/ 16389371c9d4SSatish Balay PetscErrorCode PCASMGetLocalSubmatrices(PC pc, PetscInt *n, Mat *mat[]) { 16394b9ad928SBarry Smith PC_ASM *osm; 1640ace3abfcSBarry Smith PetscBool match; 16414b9ad928SBarry Smith 16424b9ad928SBarry Smith PetscFunctionBegin; 16430700a824SBarry Smith PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1644f17a6784SPierre Jolivet if (n) PetscValidIntPointer(n, 2); 164592bb6962SLisandro Dalcin if (mat) PetscValidPointer(mat, 3); 164628b400f6SJacob Faibussowitsch PetscCheck(pc->setupcalled, PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONGSTATE, "Must call after KSPSetUp() or PCSetUp()."); 16479566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCASM, &match)); 164892bb6962SLisandro Dalcin if (!match) { 164992bb6962SLisandro Dalcin if (n) *n = 0; 16500298fd71SBarry Smith if (mat) *mat = NULL; 165192bb6962SLisandro Dalcin } else { 16524b9ad928SBarry Smith osm = (PC_ASM *)pc->data; 16534b9ad928SBarry Smith if (n) *n = osm->n_local_true; 16544b9ad928SBarry Smith if (mat) *mat = osm->pmat; 165592bb6962SLisandro Dalcin } 16564b9ad928SBarry Smith PetscFunctionReturn(0); 16574b9ad928SBarry Smith } 1658d709ea83SDmitry Karpeev 1659d709ea83SDmitry Karpeev /*@ 1660*f1580f4eSBarry Smith PCASMSetDMSubdomains - Indicates whether to use `DMCreateDomainDecomposition()` to define the subdomains, whenever possible. 1661f1ee410cSBarry Smith 1662d709ea83SDmitry Karpeev Logically Collective 1663d709ea83SDmitry Karpeev 1664d8d19677SJose E. Roman Input Parameters: 1665d709ea83SDmitry Karpeev + pc - the preconditioner 1666*f1580f4eSBarry Smith - flg - boolean indicating whether to use subdomains defined by the `DM` 1667d709ea83SDmitry Karpeev 1668d709ea83SDmitry Karpeev Options Database Key: 1669*f1580f4eSBarry Smith . -pc_asm_dm_subdomains <bool> - use subdomains defined by the `DM` 1670d709ea83SDmitry Karpeev 1671d709ea83SDmitry Karpeev Level: intermediate 1672d709ea83SDmitry Karpeev 1673*f1580f4eSBarry Smith Note: 1674*f1580f4eSBarry Smith `PCASMSetTotalSubdomains()` and `PCASMSetOverlap()` take precedence over `PCASMSetDMSubdomains()`, 1675d709ea83SDmitry Karpeev so setting either of the first two effectively turns the latter off. 1676d709ea83SDmitry Karpeev 1677*f1580f4eSBarry Smith .seealso: `PCASM`, `PCASMGetDMSubdomains()`, `PCASMSetTotalSubdomains()`, `PCASMSetOverlap()` 1678db781477SPatrick Sanan `PCASMCreateSubdomains2D()`, `PCASMSetLocalSubdomains()`, `PCASMGetLocalSubdomains()` 1679d709ea83SDmitry Karpeev @*/ 16809371c9d4SSatish Balay PetscErrorCode PCASMSetDMSubdomains(PC pc, PetscBool flg) { 1681d709ea83SDmitry Karpeev PC_ASM *osm = (PC_ASM *)pc->data; 1682d709ea83SDmitry Karpeev PetscBool match; 1683d709ea83SDmitry Karpeev 1684d709ea83SDmitry Karpeev PetscFunctionBegin; 1685d709ea83SDmitry Karpeev PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1686d709ea83SDmitry Karpeev PetscValidLogicalCollectiveBool(pc, flg, 2); 168728b400f6SJacob Faibussowitsch PetscCheck(!pc->setupcalled, ((PetscObject)pc)->comm, PETSC_ERR_ARG_WRONGSTATE, "Not for a setup PC."); 16889566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCASM, &match)); 1689ad540459SPierre Jolivet if (match) osm->dm_subdomains = flg; 1690d709ea83SDmitry Karpeev PetscFunctionReturn(0); 1691d709ea83SDmitry Karpeev } 1692d709ea83SDmitry Karpeev 1693d709ea83SDmitry Karpeev /*@ 1694*f1580f4eSBarry Smith PCASMGetDMSubdomains - Returns flag indicating whether to use `DMCreateDomainDecomposition()` to define the subdomains, whenever possible. 1695*f1580f4eSBarry Smith 1696d709ea83SDmitry Karpeev Not Collective 1697d709ea83SDmitry Karpeev 1698d709ea83SDmitry Karpeev Input Parameter: 1699d709ea83SDmitry Karpeev . pc - the preconditioner 1700d709ea83SDmitry Karpeev 1701d709ea83SDmitry Karpeev Output Parameter: 1702d709ea83SDmitry Karpeev . flg - boolean indicating whether to use subdomains defined by the DM 1703d709ea83SDmitry Karpeev 1704d709ea83SDmitry Karpeev Level: intermediate 1705d709ea83SDmitry Karpeev 1706*f1580f4eSBarry Smith .seealso: `PCASM`, `PCASMSetDMSubdomains()`, `PCASMSetTotalSubdomains()`, `PCASMSetOverlap()` 1707db781477SPatrick Sanan `PCASMCreateSubdomains2D()`, `PCASMSetLocalSubdomains()`, `PCASMGetLocalSubdomains()` 1708d709ea83SDmitry Karpeev @*/ 17099371c9d4SSatish Balay PetscErrorCode PCASMGetDMSubdomains(PC pc, PetscBool *flg) { 1710d709ea83SDmitry Karpeev PC_ASM *osm = (PC_ASM *)pc->data; 1711d709ea83SDmitry Karpeev PetscBool match; 1712d709ea83SDmitry Karpeev 1713d709ea83SDmitry Karpeev PetscFunctionBegin; 1714d709ea83SDmitry Karpeev PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1715534a8f05SLisandro Dalcin PetscValidBoolPointer(flg, 2); 17169566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc, PCASM, &match)); 171756ea09ceSStefano Zampini if (match) *flg = osm->dm_subdomains; 171856ea09ceSStefano Zampini else *flg = PETSC_FALSE; 1719d709ea83SDmitry Karpeev PetscFunctionReturn(0); 1720d709ea83SDmitry Karpeev } 172180ec0b7dSPatrick Sanan 172280ec0b7dSPatrick Sanan /*@ 1723*f1580f4eSBarry Smith PCASMGetSubMatType - Gets the matrix type used for `PCASM` subsolves, as a string. 172480ec0b7dSPatrick Sanan 172580ec0b7dSPatrick Sanan Not Collective 172680ec0b7dSPatrick Sanan 172780ec0b7dSPatrick Sanan Input Parameter: 1728*f1580f4eSBarry Smith . pc - the `PC` 172980ec0b7dSPatrick Sanan 173080ec0b7dSPatrick Sanan Output Parameter: 1731*f1580f4eSBarry Smith . pc_asm_sub_mat_type - name of matrix type 173280ec0b7dSPatrick Sanan 173380ec0b7dSPatrick Sanan Level: advanced 173480ec0b7dSPatrick Sanan 1735*f1580f4eSBarry Smith .seealso: `PCASM`, `PCASMSetSubMatType()`, `PCASM`, `PCSetType()`, `VecSetType()`, `MatType`, `Mat` 173680ec0b7dSPatrick Sanan @*/ 17379371c9d4SSatish Balay PetscErrorCode PCASMGetSubMatType(PC pc, MatType *sub_mat_type) { 173856ea09ceSStefano Zampini PetscFunctionBegin; 173956ea09ceSStefano Zampini PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1740cac4c232SBarry Smith PetscTryMethod(pc, "PCASMGetSubMatType_C", (PC, MatType *), (pc, sub_mat_type)); 174180ec0b7dSPatrick Sanan PetscFunctionReturn(0); 174280ec0b7dSPatrick Sanan } 174380ec0b7dSPatrick Sanan 174480ec0b7dSPatrick Sanan /*@ 1745*f1580f4eSBarry Smith PCASMSetSubMatType - Set the type of matrix used for `PCASM` subsolves 174680ec0b7dSPatrick Sanan 1747*f1580f4eSBarry Smith Collective on pc 174880ec0b7dSPatrick Sanan 174980ec0b7dSPatrick Sanan Input Parameters: 1750*f1580f4eSBarry Smith + pc - the `PC` object 1751*f1580f4eSBarry Smith - sub_mat_type - the `MatType` 175280ec0b7dSPatrick Sanan 175380ec0b7dSPatrick Sanan Options Database Key: 1754*f1580f4eSBarry Smith . -pc_asm_sub_mat_type <sub_mat_type> - Sets the matrix type used for subsolves, for example, seqaijviennacl. 1755*f1580f4eSBarry Smith If you specify a base name like aijviennacl, the corresponding sequential type is assumed. 175680ec0b7dSPatrick Sanan 1757*f1580f4eSBarry Smith Note: 175880ec0b7dSPatrick Sanan See "${PETSC_DIR}/include/petscmat.h" for available types 175980ec0b7dSPatrick Sanan 176080ec0b7dSPatrick Sanan Level: advanced 176180ec0b7dSPatrick Sanan 1762*f1580f4eSBarry Smith .seealso: `PCASM`, `PCASMGetSubMatType()`, `PCASM`, `PCSetType()`, `VecSetType()`, `MatType`, `Mat` 176380ec0b7dSPatrick Sanan @*/ 17649371c9d4SSatish Balay PetscErrorCode PCASMSetSubMatType(PC pc, MatType sub_mat_type) { 176556ea09ceSStefano Zampini PetscFunctionBegin; 176656ea09ceSStefano Zampini PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 1767cac4c232SBarry Smith PetscTryMethod(pc, "PCASMSetSubMatType_C", (PC, MatType), (pc, sub_mat_type)); 176880ec0b7dSPatrick Sanan PetscFunctionReturn(0); 176980ec0b7dSPatrick Sanan } 1770