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 166849ba73SBarry Smith static PetscErrorCode PCView_ASM(PC pc,PetscViewer viewer) 174b9ad928SBarry Smith { 1892bb6962SLisandro Dalcin PC_ASM *osm = (PC_ASM*)pc->data; 1913f74950SBarry Smith PetscMPIInt rank; 204d219a6aSLisandro Dalcin PetscInt i,bsz; 21ace3abfcSBarry Smith PetscBool iascii,isstring; 224b9ad928SBarry Smith PetscViewer sviewer; 23ed155784SPierre Jolivet PetscViewerFormat format; 24ed155784SPierre Jolivet const char *prefix; 254b9ad928SBarry Smith 264b9ad928SBarry Smith PetscFunctionBegin; 279566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii)); 289566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring)); 2932077d6dSBarry Smith if (iascii) { 303d03bb48SJed Brown char overlaps[256] = "user-defined overlap",blocks[256] = "total subdomain blocks not yet set"; 3163a3b9bcSJacob Faibussowitsch if (osm->overlap >= 0) PetscCall(PetscSNPrintf(overlaps,sizeof(overlaps),"amount of overlap = %" PetscInt_FMT,osm->overlap)); 3263a3b9bcSJacob Faibussowitsch if (osm->n > 0) PetscCall(PetscSNPrintf(blocks,sizeof(blocks),"total subdomain blocks = %" PetscInt_FMT,osm->n)); 339566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," %s, %s\n",blocks,overlaps)); 349566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," restriction/interpolation type - %s\n",PCASMTypes[osm->type])); 359566063dSJacob Faibussowitsch if (osm->dm_subdomains) PetscCall(PetscViewerASCIIPrintf(viewer," Additive Schwarz: using DM to define subdomains\n")); 369566063dSJacob Faibussowitsch if (osm->loctype != PC_COMPOSITE_ADDITIVE) PetscCall(PetscViewerASCIIPrintf(viewer," Additive Schwarz: local solve composition type - %s\n",PCCompositeTypes[osm->loctype])); 379566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank)); 389566063dSJacob Faibussowitsch PetscCall(PetscViewerGetFormat(viewer,&format)); 39ed155784SPierre Jolivet if (format != PETSC_VIEWER_ASCII_INFO_DETAIL) { 404d219a6aSLisandro Dalcin if (osm->ksp) { 419566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Local solver information for first block is in the following KSP and PC objects on rank 0:\n")); 429566063dSJacob Faibussowitsch PetscCall(PCGetOptionsPrefix(pc,&prefix)); 439566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Use -%sksp_view ::ascii_info_detail to display information for all blocks\n",prefix?prefix:"")); 449566063dSJacob Faibussowitsch PetscCall(PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer)); 45dd400576SPatrick Sanan if (rank == 0) { 469566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 479566063dSJacob Faibussowitsch PetscCall(KSPView(osm->ksp[0],sviewer)); 489566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 494b9ad928SBarry Smith } 509566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer)); 514d219a6aSLisandro Dalcin } 524b9ad928SBarry Smith } else { 539566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 5463a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(viewer," [%d] number of local blocks = %" PetscInt_FMT "\n",(int)rank,osm->n_local_true)); 559566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 569566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Local solver information for each block is in the following KSP and PC objects:\n")); 579566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 589566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"- - - - - - - - - - - - - - - - - -\n")); 599566063dSJacob Faibussowitsch PetscCall(PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer)); 60dd2fa690SBarry Smith for (i=0; i<osm->n_local_true; i++) { 619566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(osm->is[i],&bsz)); 6263a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(sviewer,"[%d] local block number %" PetscInt_FMT ", size = %" PetscInt_FMT "\n",(int)rank,i,bsz)); 639566063dSJacob Faibussowitsch PetscCall(KSPView(osm->ksp[i],sviewer)); 649566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(sviewer,"- - - - - - - - - - - - - - - - - -\n")); 654b9ad928SBarry Smith } 669566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer)); 679566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 689566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 699566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 704b9ad928SBarry Smith } 714b9ad928SBarry Smith } else if (isstring) { 7263a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerStringSPrintf(viewer," blocks=%" PetscInt_FMT ", overlap=%" PetscInt_FMT ", type=%s",osm->n,osm->overlap,PCASMTypes[osm->type])); 739566063dSJacob Faibussowitsch PetscCall(PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer)); 749566063dSJacob Faibussowitsch if (osm->ksp) PetscCall(KSPView(osm->ksp[0],sviewer)); 759566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer)); 764b9ad928SBarry Smith } 774b9ad928SBarry Smith PetscFunctionReturn(0); 784b9ad928SBarry Smith } 794b9ad928SBarry Smith 8092bb6962SLisandro Dalcin static PetscErrorCode PCASMPrintSubdomains(PC pc) 8192bb6962SLisandro Dalcin { 8292bb6962SLisandro Dalcin PC_ASM *osm = (PC_ASM*)pc->data; 8392bb6962SLisandro Dalcin const char *prefix; 8492bb6962SLisandro Dalcin char fname[PETSC_MAX_PATH_LEN+1]; 85643535aeSDmitry Karpeev PetscViewer viewer, sviewer; 86643535aeSDmitry Karpeev char *s; 8792bb6962SLisandro Dalcin PetscInt i,j,nidx; 8892bb6962SLisandro Dalcin const PetscInt *idx; 89643535aeSDmitry Karpeev PetscMPIInt rank, size; 90846783a0SBarry Smith 9192bb6962SLisandro Dalcin PetscFunctionBegin; 929566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size)); 939566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank)); 949566063dSJacob Faibussowitsch PetscCall(PCGetOptionsPrefix(pc,&prefix)); 959566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL,prefix,"-pc_asm_print_subdomains",fname,sizeof(fname),NULL)); 969566063dSJacob Faibussowitsch if (fname[0] == 0) PetscCall(PetscStrcpy(fname,"stdout")); 979566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIOpen(PetscObjectComm((PetscObject)pc),fname,&viewer)); 98643535aeSDmitry Karpeev for (i=0; i<osm->n_local; i++) { 99643535aeSDmitry Karpeev if (i < osm->n_local_true) { 1009566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(osm->is[i],&nidx)); 1019566063dSJacob Faibussowitsch PetscCall(ISGetIndices(osm->is[i],&idx)); 102643535aeSDmitry Karpeev /* Print to a string viewer; no more than 15 characters per index plus 512 char for the header.*/ 10336a9e3b9SBarry Smith #define len 16*(nidx+1)+512 1049566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(len,&s)); 1059566063dSJacob Faibussowitsch PetscCall(PetscViewerStringOpen(PETSC_COMM_SELF, s, len, &sviewer)); 10636a9e3b9SBarry Smith #undef len 10763a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerStringSPrintf(sviewer, "[%d:%d] Subdomain %" PetscInt_FMT " with overlap:\n", rank, size, i)); 10892bb6962SLisandro Dalcin for (j=0; j<nidx; j++) { 10963a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerStringSPrintf(sviewer,"%" PetscInt_FMT " ",idx[j])); 11092bb6962SLisandro Dalcin } 1119566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(osm->is[i],&idx)); 1129566063dSJacob Faibussowitsch PetscCall(PetscViewerStringSPrintf(sviewer,"\n")); 1139566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&sviewer)); 1149566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 11563a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "%s", s)); 1169566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 1179566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 1189566063dSJacob Faibussowitsch PetscCall(PetscFree(s)); 1192b691e39SMatthew Knepley if (osm->is_local) { 120643535aeSDmitry Karpeev /* Print to a string viewer; no more than 15 characters per index plus 512 char for the header.*/ 12136a9e3b9SBarry Smith #define len 16*(nidx+1)+512 1229566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(len, &s)); 1239566063dSJacob Faibussowitsch PetscCall(PetscViewerStringOpen(PETSC_COMM_SELF, s, len, &sviewer)); 12436a9e3b9SBarry Smith #undef len 12563a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerStringSPrintf(sviewer, "[%d:%d] Subdomain %" PetscInt_FMT " without overlap:\n", rank, size, i)); 1269566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(osm->is_local[i],&nidx)); 1279566063dSJacob Faibussowitsch PetscCall(ISGetIndices(osm->is_local[i],&idx)); 1282b691e39SMatthew Knepley for (j=0; j<nidx; j++) { 12963a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerStringSPrintf(sviewer,"%" PetscInt_FMT " ",idx[j])); 1302b691e39SMatthew Knepley } 1319566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(osm->is_local[i],&idx)); 1329566063dSJacob Faibussowitsch PetscCall(PetscViewerStringSPrintf(sviewer,"\n")); 1339566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&sviewer)); 1349566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 13563a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "%s", s)); 1369566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 1379566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 1389566063dSJacob Faibussowitsch PetscCall(PetscFree(s)); 139643535aeSDmitry Karpeev } 1402fa5cd67SKarl Rupp } else { 141643535aeSDmitry Karpeev /* Participate in collective viewer calls. */ 1429566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 1439566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 1449566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 145643535aeSDmitry Karpeev /* Assume either all ranks have is_local or none do. */ 146643535aeSDmitry Karpeev if (osm->is_local) { 1479566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 1489566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 1499566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 150643535aeSDmitry Karpeev } 151fdc48646SDmitry Karpeev } 15292bb6962SLisandro Dalcin } 1539566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 1549566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 15592bb6962SLisandro Dalcin PetscFunctionReturn(0); 15692bb6962SLisandro Dalcin } 15792bb6962SLisandro Dalcin 1586849ba73SBarry Smith static PetscErrorCode PCSetUp_ASM(PC pc) 1594b9ad928SBarry Smith { 1604b9ad928SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 16157501b6eSBarry Smith PetscBool flg; 16287e86712SBarry Smith PetscInt i,m,m_local; 1634b9ad928SBarry Smith MatReuse scall = MAT_REUSE_MATRIX; 1644b9ad928SBarry Smith IS isl; 1654b9ad928SBarry Smith KSP ksp; 1664b9ad928SBarry Smith PC subpc; 1672dcb1b2aSMatthew Knepley const char *prefix,*pprefix; 16823ce1328SBarry Smith Vec vec; 1690298fd71SBarry Smith DM *domain_dm = NULL; 1704b9ad928SBarry Smith 1714b9ad928SBarry Smith PetscFunctionBegin; 1724b9ad928SBarry Smith if (!pc->setupcalled) { 173265a2120SBarry Smith PetscInt m; 17492bb6962SLisandro Dalcin 1752b837212SDmitry Karpeev /* Note: if subdomains have been set either via PCASMSetTotalSubdomains() or via PCASMSetLocalSubdomains(), osm->n_local_true will not be PETSC_DECIDE */ 1762b837212SDmitry Karpeev if (osm->n_local_true == PETSC_DECIDE) { 17769ca1f37SDmitry Karpeev /* no subdomains given */ 17865db9045SDmitry Karpeev /* try pc->dm first, if allowed */ 179d709ea83SDmitry Karpeev if (osm->dm_subdomains && pc->dm) { 180feb221f8SDmitry Karpeev PetscInt num_domains, d; 181feb221f8SDmitry Karpeev char **domain_names; 1828d4ac253SDmitry Karpeev IS *inner_domain_is, *outer_domain_is; 1839566063dSJacob Faibussowitsch PetscCall(DMCreateDomainDecomposition(pc->dm, &num_domains, &domain_names, &inner_domain_is, &outer_domain_is, &domain_dm)); 184704f0395SPatrick Sanan osm->overlap = -1; /* We do not want to increase the overlap of the IS. 185704f0395SPatrick Sanan A future improvement of this code might allow one to use 186704f0395SPatrick Sanan DM-defined subdomains and also increase the overlap, 187704f0395SPatrick Sanan but that is not currently supported */ 188*1baa6e33SBarry Smith if (num_domains) PetscCall(PCASMSetLocalSubdomains(pc, num_domains, outer_domain_is, inner_domain_is)); 189feb221f8SDmitry Karpeev for (d = 0; d < num_domains; ++d) { 1909566063dSJacob Faibussowitsch if (domain_names) PetscCall(PetscFree(domain_names[d])); 1919566063dSJacob Faibussowitsch if (inner_domain_is) PetscCall(ISDestroy(&inner_domain_is[d])); 1929566063dSJacob Faibussowitsch if (outer_domain_is) PetscCall(ISDestroy(&outer_domain_is[d])); 193feb221f8SDmitry Karpeev } 1949566063dSJacob Faibussowitsch PetscCall(PetscFree(domain_names)); 1959566063dSJacob Faibussowitsch PetscCall(PetscFree(inner_domain_is)); 1969566063dSJacob Faibussowitsch PetscCall(PetscFree(outer_domain_is)); 197feb221f8SDmitry Karpeev } 1982b837212SDmitry Karpeev if (osm->n_local_true == PETSC_DECIDE) { 19969ca1f37SDmitry Karpeev /* still no subdomains; use one subdomain per processor */ 2002b837212SDmitry Karpeev osm->n_local_true = 1; 201feb221f8SDmitry Karpeev } 2022b837212SDmitry Karpeev } 2032b837212SDmitry Karpeev { /* determine the global and max number of subdomains */ 2046ac3741eSJed Brown struct {PetscInt max,sum;} inwork,outwork; 205c10200c1SHong Zhang PetscMPIInt size; 206c10200c1SHong Zhang 2076ac3741eSJed Brown inwork.max = osm->n_local_true; 2086ac3741eSJed Brown inwork.sum = osm->n_local_true; 2091c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(&inwork,&outwork,1,MPIU_2INT,MPIU_MAXSUM_OP,PetscObjectComm((PetscObject)pc))); 2106ac3741eSJed Brown osm->n_local = outwork.max; 2116ac3741eSJed Brown osm->n = outwork.sum; 212c10200c1SHong Zhang 2139566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size)); 214c10200c1SHong Zhang if (outwork.max == 1 && outwork.sum == size) { 2157dae84e0SHong Zhang /* osm->n_local_true = 1 on all processes, set this option may enable use of optimized MatCreateSubMatrices() implementation */ 2169566063dSJacob Faibussowitsch PetscCall(MatSetOption(pc->pmat,MAT_SUBMAT_SINGLEIS,PETSC_TRUE)); 217c10200c1SHong Zhang } 2184b9ad928SBarry Smith } 21978904715SLisandro Dalcin if (!osm->is) { /* create the index sets */ 2209566063dSJacob Faibussowitsch PetscCall(PCASMCreateSubdomains(pc->pmat,osm->n_local_true,&osm->is)); 2214b9ad928SBarry Smith } 222f5234e35SJed Brown if (osm->n_local_true > 1 && !osm->is_local) { 2239566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(osm->n_local_true,&osm->is_local)); 224f5234e35SJed Brown for (i=0; i<osm->n_local_true; i++) { 225f5234e35SJed Brown if (osm->overlap > 0) { /* With positive overlap, osm->is[i] will be modified */ 2269566063dSJacob Faibussowitsch PetscCall(ISDuplicate(osm->is[i],&osm->is_local[i])); 2279566063dSJacob Faibussowitsch PetscCall(ISCopy(osm->is[i],osm->is_local[i])); 228f5234e35SJed Brown } else { 2299566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)osm->is[i])); 230f5234e35SJed Brown osm->is_local[i] = osm->is[i]; 231f5234e35SJed Brown } 232f5234e35SJed Brown } 233f5234e35SJed Brown } 2349566063dSJacob Faibussowitsch PetscCall(PCGetOptionsPrefix(pc,&prefix)); 2353d03bb48SJed Brown if (osm->overlap > 0) { 2364b9ad928SBarry Smith /* Extend the "overlapping" regions by a number of steps */ 2379566063dSJacob Faibussowitsch PetscCall(MatIncreaseOverlap(pc->pmat,osm->n_local_true,osm->is,osm->overlap)); 2383d03bb48SJed Brown } 2396ed231c7SMatthew Knepley if (osm->sort_indices) { 24092bb6962SLisandro Dalcin for (i=0; i<osm->n_local_true; i++) { 2419566063dSJacob Faibussowitsch PetscCall(ISSort(osm->is[i])); 2422b691e39SMatthew Knepley if (osm->is_local) { 2439566063dSJacob Faibussowitsch PetscCall(ISSort(osm->is_local[i])); 2442b691e39SMatthew Knepley } 2454b9ad928SBarry Smith } 2466ed231c7SMatthew Knepley } 24798d3e228SPierre Jolivet flg = PETSC_FALSE; 24898d3e228SPierre Jolivet PetscCall(PetscOptionsHasName(NULL,prefix,"-pc_asm_print_subdomains",&flg)); 24998d3e228SPierre Jolivet if (flg) PetscCall(PCASMPrintSubdomains(pc)); 250f6991133SBarry Smith if (!osm->ksp) { 25178904715SLisandro Dalcin /* Create the local solvers */ 2529566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(osm->n_local_true,&osm->ksp)); 253feb221f8SDmitry Karpeev if (domain_dm) { 2549566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc,"Setting up ASM subproblems using the embedded DM\n")); 255feb221f8SDmitry Karpeev } 25692bb6962SLisandro Dalcin for (i=0; i<osm->n_local_true; i++) { 2579566063dSJacob Faibussowitsch PetscCall(KSPCreate(PETSC_COMM_SELF,&ksp)); 2589566063dSJacob Faibussowitsch PetscCall(KSPSetErrorIfNotConverged(ksp,pc->erroriffailure)); 2599566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)ksp)); 2609566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,1)); 2619566063dSJacob Faibussowitsch PetscCall(KSPSetType(ksp,KSPPREONLY)); 2629566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp,&subpc)); 2639566063dSJacob Faibussowitsch PetscCall(PCGetOptionsPrefix(pc,&prefix)); 2649566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(ksp,prefix)); 2659566063dSJacob Faibussowitsch PetscCall(KSPAppendOptionsPrefix(ksp,"sub_")); 266feb221f8SDmitry Karpeev if (domain_dm) { 2679566063dSJacob Faibussowitsch PetscCall(KSPSetDM(ksp, domain_dm[i])); 2689566063dSJacob Faibussowitsch PetscCall(KSPSetDMActive(ksp, PETSC_FALSE)); 2699566063dSJacob Faibussowitsch PetscCall(DMDestroy(&domain_dm[i])); 270feb221f8SDmitry Karpeev } 2714b9ad928SBarry Smith osm->ksp[i] = ksp; 2724b9ad928SBarry Smith } 273*1baa6e33SBarry Smith if (domain_dm) PetscCall(PetscFree(domain_dm)); 274f6991133SBarry Smith } 2751dd8081eSeaulisa 2769566063dSJacob Faibussowitsch PetscCall(ISConcatenate(PETSC_COMM_SELF, osm->n_local_true, osm->is, &osm->lis)); 2779566063dSJacob Faibussowitsch PetscCall(ISSortRemoveDups(osm->lis)); 2789566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(osm->lis, &m)); 2791dd8081eSeaulisa 2804b9ad928SBarry Smith scall = MAT_INITIAL_MATRIX; 2814b9ad928SBarry Smith } else { 2824b9ad928SBarry Smith /* 2834b9ad928SBarry Smith Destroy the blocks from the previous iteration 2844b9ad928SBarry Smith */ 285e09e08ccSBarry Smith if (pc->flag == DIFFERENT_NONZERO_PATTERN) { 2869566063dSJacob Faibussowitsch PetscCall(MatDestroyMatrices(osm->n_local_true,&osm->pmat)); 2874b9ad928SBarry Smith scall = MAT_INITIAL_MATRIX; 2884b9ad928SBarry Smith } 2894b9ad928SBarry Smith } 2904b9ad928SBarry Smith 291b58cb649SBarry Smith /* Destroy previous submatrices of a different type than pc->pmat since MAT_REUSE_MATRIX won't work in that case */ 29282b5ce2aSStefano Zampini if (scall == MAT_REUSE_MATRIX && osm->sub_mat_type) { 293b58cb649SBarry Smith if (osm->n_local_true > 0) { 2949566063dSJacob Faibussowitsch PetscCall(MatDestroySubMatrices(osm->n_local_true,&osm->pmat)); 295b58cb649SBarry Smith } 296b58cb649SBarry Smith scall = MAT_INITIAL_MATRIX; 297b58cb649SBarry Smith } 298b58cb649SBarry Smith 29992bb6962SLisandro Dalcin /* 30092bb6962SLisandro Dalcin Extract out the submatrices 30192bb6962SLisandro Dalcin */ 3029566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(pc->pmat,osm->n_local_true,osm->is,osm->is,scall,&osm->pmat)); 30392bb6962SLisandro Dalcin if (scall == MAT_INITIAL_MATRIX) { 3049566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptionsPrefix((PetscObject)pc->pmat,&pprefix)); 30592bb6962SLisandro Dalcin for (i=0; i<osm->n_local_true; i++) { 3069566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)osm->pmat[i])); 3079566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)osm->pmat[i],pprefix)); 30892bb6962SLisandro Dalcin } 30992bb6962SLisandro Dalcin } 31080ec0b7dSPatrick Sanan 31180ec0b7dSPatrick Sanan /* Convert the types of the submatrices (if needbe) */ 31280ec0b7dSPatrick Sanan if (osm->sub_mat_type) { 31380ec0b7dSPatrick Sanan for (i=0; i<osm->n_local_true; i++) { 3149566063dSJacob Faibussowitsch PetscCall(MatConvert(osm->pmat[i],osm->sub_mat_type,MAT_INPLACE_MATRIX,&(osm->pmat[i]))); 31580ec0b7dSPatrick Sanan } 31680ec0b7dSPatrick Sanan } 31780ec0b7dSPatrick Sanan 31880ec0b7dSPatrick Sanan if (!pc->setupcalled) { 31956ea09ceSStefano Zampini VecType vtype; 32056ea09ceSStefano Zampini 32180ec0b7dSPatrick Sanan /* Create the local work vectors (from the local matrices) and scatter contexts */ 3229566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(pc->pmat,&vec,NULL)); 3235b3c0d42Seaulisa 3247827d75bSBarry 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()"); 3251dd8081eSeaulisa if (osm->is_local && osm->type == PC_ASM_RESTRICT && osm->loctype == PC_COMPOSITE_ADDITIVE) { 3269566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(osm->n_local_true,&osm->lprolongation)); 3271dd8081eSeaulisa } 3289566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(osm->n_local_true,&osm->lrestriction)); 3299566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(osm->n_local_true,&osm->x)); 3309566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(osm->n_local_true,&osm->y)); 331347574c9Seaulisa 3329566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(osm->lis,&m)); 3339566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF,m,0,1,&isl)); 3349566063dSJacob Faibussowitsch PetscCall(MatGetVecType(osm->pmat[0],&vtype)); 3359566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_SELF,&osm->lx)); 3369566063dSJacob Faibussowitsch PetscCall(VecSetSizes(osm->lx,m,m)); 3379566063dSJacob Faibussowitsch PetscCall(VecSetType(osm->lx,vtype)); 3389566063dSJacob Faibussowitsch PetscCall(VecDuplicate(osm->lx, &osm->ly)); 3399566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(vec,osm->lis,osm->lx,isl,&osm->restriction)); 3409566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isl)); 341347574c9Seaulisa 34280ec0b7dSPatrick Sanan for (i=0; i<osm->n_local_true; ++i) { 3435b3c0d42Seaulisa ISLocalToGlobalMapping ltog; 3445b3c0d42Seaulisa IS isll; 3455b3c0d42Seaulisa const PetscInt *idx_is; 3465b3c0d42Seaulisa PetscInt *idx_lis,nout; 3475b3c0d42Seaulisa 3489566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(osm->is[i],&m)); 3499566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(osm->pmat[i],&osm->x[i],NULL)); 3509566063dSJacob Faibussowitsch PetscCall(VecDuplicate(osm->x[i],&osm->y[i])); 35155817e79SBarry Smith 352b0de9f23SBarry Smith /* generate a scatter from ly to y[i] picking all the overlapping is[i] entries */ 3539566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreateIS(osm->lis,<og)); 3549566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(osm->is[i],&m)); 3559566063dSJacob Faibussowitsch PetscCall(ISGetIndices(osm->is[i], &idx_is)); 3569566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m,&idx_lis)); 3579566063dSJacob Faibussowitsch PetscCall(ISGlobalToLocalMappingApply(ltog,IS_GTOLM_DROP,m,idx_is,&nout,idx_lis)); 3587827d75bSBarry Smith PetscCheck(nout == m,PETSC_COMM_SELF,PETSC_ERR_PLIB,"is not a subset of lis"); 3599566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(osm->is[i], &idx_is)); 3609566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF,m,idx_lis,PETSC_OWN_POINTER,&isll)); 3619566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingDestroy(<og)); 3629566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF,m,0,1,&isl)); 3639566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(osm->ly,isll,osm->y[i],isl,&osm->lrestriction[i])); 3649566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isll)); 3659566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isl)); 366910cf402Sprj- if (osm->lprolongation) { /* generate a scatter from y[i] to ly picking only the the non-overlapping is_local[i] entries */ 367d8b3b5e3Seaulisa ISLocalToGlobalMapping ltog; 368d8b3b5e3Seaulisa IS isll,isll_local; 369d8b3b5e3Seaulisa const PetscInt *idx_local; 370d8b3b5e3Seaulisa PetscInt *idx1, *idx2, nout; 371d8b3b5e3Seaulisa 3729566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(osm->is_local[i],&m_local)); 3739566063dSJacob Faibussowitsch PetscCall(ISGetIndices(osm->is_local[i], &idx_local)); 374d8b3b5e3Seaulisa 3759566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreateIS(osm->is[i],<og)); 3769566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m_local,&idx1)); 3779566063dSJacob Faibussowitsch PetscCall(ISGlobalToLocalMappingApply(ltog,IS_GTOLM_DROP,m_local,idx_local,&nout,idx1)); 3789566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingDestroy(<og)); 3797827d75bSBarry Smith PetscCheck(nout == m_local,PETSC_COMM_SELF,PETSC_ERR_PLIB,"is_local not a subset of is"); 3809566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF,m_local,idx1,PETSC_OWN_POINTER,&isll)); 381d8b3b5e3Seaulisa 3829566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreateIS(osm->lis,<og)); 3839566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m_local,&idx2)); 3849566063dSJacob Faibussowitsch PetscCall(ISGlobalToLocalMappingApply(ltog,IS_GTOLM_DROP,m_local,idx_local,&nout,idx2)); 3859566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingDestroy(<og)); 3867827d75bSBarry Smith PetscCheck(nout == m_local,PETSC_COMM_SELF,PETSC_ERR_PLIB,"is_local not a subset of lis"); 3879566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF,m_local,idx2,PETSC_OWN_POINTER,&isll_local)); 388d8b3b5e3Seaulisa 3899566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(osm->is_local[i], &idx_local)); 3909566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(osm->y[i],isll,osm->ly,isll_local,&osm->lprolongation[i])); 391d8b3b5e3Seaulisa 3929566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isll)); 3939566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isll_local)); 394d8b3b5e3Seaulisa } 39580ec0b7dSPatrick Sanan } 3969566063dSJacob Faibussowitsch PetscCall(VecDestroy(&vec)); 39780ec0b7dSPatrick Sanan } 39880ec0b7dSPatrick Sanan 399fb745f2cSMatthew G. Knepley if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE) { 400235cc4ceSMatthew G. Knepley IS *cis; 401235cc4ceSMatthew G. Knepley PetscInt c; 402235cc4ceSMatthew G. Knepley 4039566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(osm->n_local_true, &cis)); 404235cc4ceSMatthew G. Knepley for (c = 0; c < osm->n_local_true; ++c) cis[c] = osm->lis; 4059566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(pc->pmat, osm->n_local_true, osm->is, cis, scall, &osm->lmats)); 4069566063dSJacob Faibussowitsch PetscCall(PetscFree(cis)); 407fb745f2cSMatthew G. Knepley } 4084b9ad928SBarry Smith 4094b9ad928SBarry Smith /* Return control to the user so that the submatrices can be modified (e.g., to apply 4104b9ad928SBarry Smith different boundary conditions for the submatrices than for the global problem) */ 4119566063dSJacob Faibussowitsch PetscCall(PCModifySubMatrices(pc,osm->n_local_true,osm->is,osm->is,osm->pmat,pc->modifysubmatricesP)); 4124b9ad928SBarry Smith 41392bb6962SLisandro Dalcin /* 41492bb6962SLisandro Dalcin Loop over subdomains putting them into local ksp 41592bb6962SLisandro Dalcin */ 4169566063dSJacob Faibussowitsch PetscCall(KSPGetOptionsPrefix(osm->ksp[0],&prefix)); 41792bb6962SLisandro Dalcin for (i=0; i<osm->n_local_true; i++) { 4189566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(osm->ksp[i],osm->pmat[i],osm->pmat[i])); 4199566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(osm->pmat[i],prefix)); 42092bb6962SLisandro Dalcin if (!pc->setupcalled) { 4219566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(osm->ksp[i])); 4224b9ad928SBarry Smith } 42392bb6962SLisandro Dalcin } 4244b9ad928SBarry Smith PetscFunctionReturn(0); 4254b9ad928SBarry Smith } 4264b9ad928SBarry Smith 4276849ba73SBarry Smith static PetscErrorCode PCSetUpOnBlocks_ASM(PC pc) 4284b9ad928SBarry Smith { 4294b9ad928SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 43013f74950SBarry Smith PetscInt i; 431539c167fSBarry Smith KSPConvergedReason reason; 4324b9ad928SBarry Smith 4334b9ad928SBarry Smith PetscFunctionBegin; 4344b9ad928SBarry Smith for (i=0; i<osm->n_local_true; i++) { 4359566063dSJacob Faibussowitsch PetscCall(KSPSetUp(osm->ksp[i])); 4369566063dSJacob Faibussowitsch PetscCall(KSPGetConvergedReason(osm->ksp[i],&reason)); 437c0decd05SBarry Smith if (reason == KSP_DIVERGED_PC_FAILED) { 438261222b2SHong Zhang pc->failedreason = PC_SUBPC_ERROR; 439e0eafd54SHong Zhang } 4404b9ad928SBarry Smith } 4414b9ad928SBarry Smith PetscFunctionReturn(0); 4424b9ad928SBarry Smith } 4434b9ad928SBarry Smith 4446849ba73SBarry Smith static PetscErrorCode PCApply_ASM(PC pc,Vec x,Vec y) 4454b9ad928SBarry Smith { 4464b9ad928SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 4471dd8081eSeaulisa PetscInt i,n_local_true = osm->n_local_true; 4484b9ad928SBarry Smith ScatterMode forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE; 4494b9ad928SBarry Smith 4504b9ad928SBarry Smith PetscFunctionBegin; 4514b9ad928SBarry Smith /* 45248e38a8aSPierre Jolivet support for limiting the restriction or interpolation to only local 4534b9ad928SBarry Smith subdomain values (leaving the other values 0). 4544b9ad928SBarry Smith */ 4554b9ad928SBarry Smith if (!(osm->type & PC_ASM_RESTRICT)) { 4564b9ad928SBarry Smith forward = SCATTER_FORWARD_LOCAL; 4574b9ad928SBarry Smith /* have to zero the work RHS since scatter may leave some slots empty */ 4589566063dSJacob Faibussowitsch PetscCall(VecSet(osm->lx, 0.0)); 4594b9ad928SBarry Smith } 460347574c9Seaulisa if (!(osm->type & PC_ASM_INTERPOLATE)) { 461347574c9Seaulisa reverse = SCATTER_REVERSE_LOCAL; 462347574c9Seaulisa } 4634b9ad928SBarry Smith 464347574c9Seaulisa if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE || osm->loctype == PC_COMPOSITE_ADDITIVE) { 465b0de9f23SBarry Smith /* zero the global and the local solutions */ 4669566063dSJacob Faibussowitsch PetscCall(VecSet(y, 0.0)); 4679566063dSJacob Faibussowitsch PetscCall(VecSet(osm->ly, 0.0)); 468347574c9Seaulisa 46948e38a8aSPierre Jolivet /* copy the global RHS to local RHS including the ghost nodes */ 4709566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->restriction, x, osm->lx, INSERT_VALUES, forward)); 4719566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->restriction, x, osm->lx, INSERT_VALUES, forward)); 472347574c9Seaulisa 47348e38a8aSPierre Jolivet /* restrict local RHS to the overlapping 0-block RHS */ 4749566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward)); 4759566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward)); 476d8b3b5e3Seaulisa 47712cd4985SMatthew G. Knepley /* do the local solves */ 47812cd4985SMatthew G. Knepley for (i = 0; i < n_local_true; ++i) { 479347574c9Seaulisa 480b0de9f23SBarry Smith /* solve the overlapping i-block */ 4819566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PC_ApplyOnBlocks, osm->ksp[i], osm->x[i], osm->y[i],0)); 4829566063dSJacob Faibussowitsch PetscCall(KSPSolve(osm->ksp[i], osm->x[i], osm->y[i])); 4839566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(osm->ksp[i], pc, osm->y[i])); 4849566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PC_ApplyOnBlocks, osm->ksp[i], osm->x[i], osm->y[i], 0)); 485d8b3b5e3Seaulisa 486910cf402Sprj- if (osm->lprolongation) { /* interpolate the non-overlapping i-block solution to the local solution (only for restrictive additive) */ 4879566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward)); 4889566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward)); 48948e38a8aSPierre Jolivet } else { /* interpolate the overlapping i-block solution to the local solution */ 4909566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse)); 4919566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse)); 492d8b3b5e3Seaulisa } 493347574c9Seaulisa 494347574c9Seaulisa if (i < n_local_true-1) { 49548e38a8aSPierre Jolivet /* restrict local RHS to the overlapping (i+1)-block RHS */ 4969566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->lrestriction[i+1], osm->lx, osm->x[i+1], INSERT_VALUES, forward)); 4979566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->lrestriction[i+1], osm->lx, osm->x[i+1], INSERT_VALUES, forward)); 498347574c9Seaulisa 499347574c9Seaulisa if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE) { 5008966356dSPierre Jolivet /* update the overlapping (i+1)-block RHS using the current local solution */ 5019566063dSJacob Faibussowitsch PetscCall(MatMult(osm->lmats[i+1], osm->ly, osm->y[i+1])); 5029566063dSJacob Faibussowitsch PetscCall(VecAXPBY(osm->x[i+1],-1.,1., osm->y[i+1])); 5037c3d802fSMatthew G. Knepley } 50412cd4985SMatthew G. Knepley } 50512cd4985SMatthew G. Knepley } 50648e38a8aSPierre Jolivet /* add the local solution to the global solution including the ghost nodes */ 5079566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->restriction, osm->ly, y, ADD_VALUES, reverse)); 5089566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->restriction, osm->ly, y, ADD_VALUES, reverse)); 50998921bdaSJacob Faibussowitsch } else SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Invalid local composition type: %s", PCCompositeTypes[osm->loctype]); 5104b9ad928SBarry Smith PetscFunctionReturn(0); 5114b9ad928SBarry Smith } 5124b9ad928SBarry Smith 51348e38a8aSPierre Jolivet static PetscErrorCode PCMatApply_ASM(PC pc,Mat X,Mat Y) 51448e38a8aSPierre Jolivet { 51548e38a8aSPierre Jolivet PC_ASM *osm = (PC_ASM*)pc->data; 51648e38a8aSPierre Jolivet Mat Z,W; 51748e38a8aSPierre Jolivet Vec x; 51848e38a8aSPierre Jolivet PetscInt i,m,N; 51948e38a8aSPierre Jolivet ScatterMode forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE; 52048e38a8aSPierre Jolivet 52148e38a8aSPierre Jolivet PetscFunctionBegin; 5227827d75bSBarry Smith PetscCheck(osm->n_local_true <= 1,PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Not yet implemented"); 52348e38a8aSPierre Jolivet /* 52448e38a8aSPierre Jolivet support for limiting the restriction or interpolation to only local 52548e38a8aSPierre Jolivet subdomain values (leaving the other values 0). 52648e38a8aSPierre Jolivet */ 52748e38a8aSPierre Jolivet if (!(osm->type & PC_ASM_RESTRICT)) { 52848e38a8aSPierre Jolivet forward = SCATTER_FORWARD_LOCAL; 52948e38a8aSPierre Jolivet /* have to zero the work RHS since scatter may leave some slots empty */ 5309566063dSJacob Faibussowitsch PetscCall(VecSet(osm->lx, 0.0)); 53148e38a8aSPierre Jolivet } 53248e38a8aSPierre Jolivet if (!(osm->type & PC_ASM_INTERPOLATE)) { 53348e38a8aSPierre Jolivet reverse = SCATTER_REVERSE_LOCAL; 53448e38a8aSPierre Jolivet } 5359566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(osm->x[0], &m)); 5369566063dSJacob Faibussowitsch PetscCall(MatGetSize(X, NULL, &N)); 5379566063dSJacob Faibussowitsch PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, m, N, NULL, &Z)); 53848e38a8aSPierre Jolivet if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE || osm->loctype == PC_COMPOSITE_ADDITIVE) { 53948e38a8aSPierre Jolivet /* zero the global and the local solutions */ 5409566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(Y)); 5419566063dSJacob Faibussowitsch PetscCall(VecSet(osm->ly, 0.0)); 54248e38a8aSPierre Jolivet 54348e38a8aSPierre Jolivet for (i = 0; i < N; ++i) { 5449566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumnVecRead(X, i, &x)); 54548e38a8aSPierre Jolivet /* copy the global RHS to local RHS including the ghost nodes */ 5469566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->restriction, x, osm->lx, INSERT_VALUES, forward)); 5479566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->restriction, x, osm->lx, INSERT_VALUES, forward)); 5489566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumnVecRead(X, i, &x)); 54948e38a8aSPierre Jolivet 5509566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumnVecWrite(Z, i, &x)); 55148e38a8aSPierre Jolivet /* restrict local RHS to the overlapping 0-block RHS */ 5529566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->lrestriction[0], osm->lx, x, INSERT_VALUES, forward)); 5539566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->lrestriction[0], osm->lx, x, INSERT_VALUES, forward)); 5549566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumnVecWrite(Z, i, &x)); 55548e38a8aSPierre Jolivet } 5569566063dSJacob Faibussowitsch PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, m, N, NULL, &W)); 55748e38a8aSPierre Jolivet /* solve the overlapping 0-block */ 5589566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PC_ApplyOnBlocks, osm->ksp[0], Z, W, 0)); 5599566063dSJacob Faibussowitsch PetscCall(KSPMatSolve(osm->ksp[0], Z, W)); 5609566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(osm->ksp[0], pc, NULL)); 5619566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PC_ApplyOnBlocks, osm->ksp[0], Z, W,0)); 5629566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Z)); 56348e38a8aSPierre Jolivet 56448e38a8aSPierre Jolivet for (i = 0; i < N; ++i) { 5659566063dSJacob Faibussowitsch PetscCall(VecSet(osm->ly, 0.0)); 5669566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumnVecRead(W, i, &x)); 56748e38a8aSPierre Jolivet if (osm->lprolongation) { /* interpolate the non-overlapping 0-block solution to the local solution (only for restrictive additive) */ 5689566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->lprolongation[0], x, osm->ly, ADD_VALUES, forward)); 5699566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->lprolongation[0], x, osm->ly, ADD_VALUES, forward)); 57048e38a8aSPierre Jolivet } else { /* interpolate the overlapping 0-block solution to the local solution */ 5719566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->lrestriction[0], x, osm->ly, ADD_VALUES, reverse)); 5729566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->lrestriction[0], x, osm->ly, ADD_VALUES, reverse)); 57348e38a8aSPierre Jolivet } 5749566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumnVecRead(W, i, &x)); 57548e38a8aSPierre Jolivet 5769566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumnVecWrite(Y, i, &x)); 57748e38a8aSPierre Jolivet /* add the local solution to the global solution including the ghost nodes */ 5789566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->restriction, osm->ly, x, ADD_VALUES, reverse)); 5799566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->restriction, osm->ly, x, ADD_VALUES, reverse)); 5809566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumnVecWrite(Y, i, &x)); 58148e38a8aSPierre Jolivet } 5829566063dSJacob Faibussowitsch PetscCall(MatDestroy(&W)); 58398921bdaSJacob Faibussowitsch } else SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Invalid local composition type: %s", PCCompositeTypes[osm->loctype]); 58448e38a8aSPierre Jolivet PetscFunctionReturn(0); 58548e38a8aSPierre Jolivet } 58648e38a8aSPierre Jolivet 5876849ba73SBarry Smith static PetscErrorCode PCApplyTranspose_ASM(PC pc,Vec x,Vec y) 5884b9ad928SBarry Smith { 5894b9ad928SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 5901dd8081eSeaulisa PetscInt i,n_local_true = osm->n_local_true; 5914b9ad928SBarry Smith ScatterMode forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE; 5924b9ad928SBarry Smith 5934b9ad928SBarry Smith PetscFunctionBegin; 5944b9ad928SBarry Smith /* 5954b9ad928SBarry Smith Support for limiting the restriction or interpolation to only local 5964b9ad928SBarry Smith subdomain values (leaving the other values 0). 5974b9ad928SBarry Smith 5984b9ad928SBarry Smith Note: these are reversed from the PCApply_ASM() because we are applying the 5994b9ad928SBarry Smith transpose of the three terms 6004b9ad928SBarry Smith */ 601d8b3b5e3Seaulisa 6024b9ad928SBarry Smith if (!(osm->type & PC_ASM_INTERPOLATE)) { 6034b9ad928SBarry Smith forward = SCATTER_FORWARD_LOCAL; 6044b9ad928SBarry Smith /* have to zero the work RHS since scatter may leave some slots empty */ 6059566063dSJacob Faibussowitsch PetscCall(VecSet(osm->lx, 0.0)); 6064b9ad928SBarry Smith } 6072fa5cd67SKarl Rupp if (!(osm->type & PC_ASM_RESTRICT)) reverse = SCATTER_REVERSE_LOCAL; 6084b9ad928SBarry Smith 609b0de9f23SBarry Smith /* zero the global and the local solutions */ 6109566063dSJacob Faibussowitsch PetscCall(VecSet(y, 0.0)); 6119566063dSJacob Faibussowitsch PetscCall(VecSet(osm->ly, 0.0)); 612d8b3b5e3Seaulisa 613b0de9f23SBarry Smith /* Copy the global RHS to local RHS including the ghost nodes */ 6149566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->restriction, x, osm->lx, INSERT_VALUES, forward)); 6159566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->restriction, x, osm->lx, INSERT_VALUES, forward)); 616d8b3b5e3Seaulisa 617b0de9f23SBarry Smith /* Restrict local RHS to the overlapping 0-block RHS */ 6189566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward)); 6199566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward)); 620d8b3b5e3Seaulisa 6214b9ad928SBarry Smith /* do the local solves */ 622d8b3b5e3Seaulisa for (i = 0; i < n_local_true; ++i) { 623d8b3b5e3Seaulisa 624b0de9f23SBarry Smith /* solve the overlapping i-block */ 6259566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PC_ApplyOnBlocks,osm->ksp[i],osm->x[i],osm->y[i],0)); 6269566063dSJacob Faibussowitsch PetscCall(KSPSolveTranspose(osm->ksp[i], osm->x[i], osm->y[i])); 6279566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(osm->ksp[i],pc,osm->y[i])); 6289566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PC_ApplyOnBlocks,osm->ksp[i],osm->x[i],osm->y[i],0)); 629d8b3b5e3Seaulisa 630910cf402Sprj- if (osm->lprolongation) { /* interpolate the non-overlapping i-block solution to the local solution */ 6319566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward)); 6329566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward)); 633910cf402Sprj- } else { /* interpolate the overlapping i-block solution to the local solution */ 6349566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse)); 6359566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse)); 6364b9ad928SBarry Smith } 637d8b3b5e3Seaulisa 638d8b3b5e3Seaulisa if (i < n_local_true-1) { 639b0de9f23SBarry Smith /* Restrict local RHS to the overlapping (i+1)-block RHS */ 6409566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->lrestriction[i+1], osm->lx, osm->x[i+1], INSERT_VALUES, forward)); 6419566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->lrestriction[i+1], osm->lx, osm->x[i+1], INSERT_VALUES, forward)); 6424b9ad928SBarry Smith } 6434b9ad928SBarry Smith } 644b0de9f23SBarry Smith /* Add the local solution to the global solution including the ghost nodes */ 6459566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->restriction, osm->ly, y, ADD_VALUES, reverse)); 6469566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->restriction, osm->ly, y, ADD_VALUES, reverse)); 6474b9ad928SBarry Smith PetscFunctionReturn(0); 6484b9ad928SBarry Smith } 6494b9ad928SBarry Smith 650e91c6855SBarry Smith static PetscErrorCode PCReset_ASM(PC pc) 6514b9ad928SBarry Smith { 6524b9ad928SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 65313f74950SBarry Smith PetscInt i; 6544b9ad928SBarry Smith 6554b9ad928SBarry Smith PetscFunctionBegin; 65692bb6962SLisandro Dalcin if (osm->ksp) { 65792bb6962SLisandro Dalcin for (i=0; i<osm->n_local_true; i++) { 6589566063dSJacob Faibussowitsch PetscCall(KSPReset(osm->ksp[i])); 65992bb6962SLisandro Dalcin } 66092bb6962SLisandro Dalcin } 661e09e08ccSBarry Smith if (osm->pmat) { 66292bb6962SLisandro Dalcin if (osm->n_local_true > 0) { 6639566063dSJacob Faibussowitsch PetscCall(MatDestroySubMatrices(osm->n_local_true,&osm->pmat)); 66492bb6962SLisandro Dalcin } 66592bb6962SLisandro Dalcin } 6661dd8081eSeaulisa if (osm->lrestriction) { 6679566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&osm->restriction)); 6681dd8081eSeaulisa for (i=0; i<osm->n_local_true; i++) { 6699566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&osm->lrestriction[i])); 6709566063dSJacob Faibussowitsch if (osm->lprolongation) PetscCall(VecScatterDestroy(&osm->lprolongation[i])); 6719566063dSJacob Faibussowitsch PetscCall(VecDestroy(&osm->x[i])); 6729566063dSJacob Faibussowitsch PetscCall(VecDestroy(&osm->y[i])); 6734b9ad928SBarry Smith } 6749566063dSJacob Faibussowitsch PetscCall(PetscFree(osm->lrestriction)); 6759566063dSJacob Faibussowitsch if (osm->lprolongation) PetscCall(PetscFree(osm->lprolongation)); 6769566063dSJacob Faibussowitsch PetscCall(PetscFree(osm->x)); 6779566063dSJacob Faibussowitsch PetscCall(PetscFree(osm->y)); 6781dd8081eSeaulisa 67992bb6962SLisandro Dalcin } 6809566063dSJacob Faibussowitsch PetscCall(PCASMDestroySubdomains(osm->n_local_true,osm->is,osm->is_local)); 6819566063dSJacob Faibussowitsch PetscCall(ISDestroy(&osm->lis)); 6829566063dSJacob Faibussowitsch PetscCall(VecDestroy(&osm->lx)); 6839566063dSJacob Faibussowitsch PetscCall(VecDestroy(&osm->ly)); 684347574c9Seaulisa if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE) { 6859566063dSJacob Faibussowitsch PetscCall(MatDestroyMatrices(osm->n_local_true, &osm->lmats)); 686fb745f2cSMatthew G. Knepley } 6872fa5cd67SKarl Rupp 6889566063dSJacob Faibussowitsch PetscCall(PetscFree(osm->sub_mat_type)); 68980ec0b7dSPatrick Sanan 6900a545947SLisandro Dalcin osm->is = NULL; 6910a545947SLisandro Dalcin osm->is_local = NULL; 692e91c6855SBarry Smith PetscFunctionReturn(0); 693e91c6855SBarry Smith } 694e91c6855SBarry Smith 695e91c6855SBarry Smith static PetscErrorCode PCDestroy_ASM(PC pc) 696e91c6855SBarry Smith { 697e91c6855SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 698e91c6855SBarry Smith PetscInt i; 699e91c6855SBarry Smith 700e91c6855SBarry Smith PetscFunctionBegin; 7019566063dSJacob Faibussowitsch PetscCall(PCReset_ASM(pc)); 702e91c6855SBarry Smith if (osm->ksp) { 703e91c6855SBarry Smith for (i=0; i<osm->n_local_true; i++) { 7049566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&osm->ksp[i])); 705e91c6855SBarry Smith } 7069566063dSJacob Faibussowitsch PetscCall(PetscFree(osm->ksp)); 707e91c6855SBarry Smith } 7089566063dSJacob Faibussowitsch PetscCall(PetscFree(pc->data)); 70996322394SPierre Jolivet 7109566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCASMSetLocalSubdomains_C",NULL)); 7119566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCASMSetTotalSubdomains_C",NULL)); 7129566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCASMSetOverlap_C",NULL)); 7139566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCASMSetType_C",NULL)); 7149566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCASMGetType_C",NULL)); 7159566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCASMSetLocalType_C",NULL)); 7169566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCASMGetLocalType_C",NULL)); 7179566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCASMSetSortIndices_C",NULL)); 7189566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCASMGetSubKSP_C",NULL)); 7199566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCASMGetSubMatType_C",NULL)); 7209566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCASMSetSubMatType_C",NULL)); 7214b9ad928SBarry Smith PetscFunctionReturn(0); 7224b9ad928SBarry Smith } 7234b9ad928SBarry Smith 7244416b707SBarry Smith static PetscErrorCode PCSetFromOptions_ASM(PetscOptionItems *PetscOptionsObject,PC pc) 7254b9ad928SBarry Smith { 7264b9ad928SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 7279dcbbd2bSBarry Smith PetscInt blocks,ovl; 72857501b6eSBarry Smith PetscBool flg; 72992bb6962SLisandro Dalcin PCASMType asmtype; 73012cd4985SMatthew G. Knepley PCCompositeType loctype; 73180ec0b7dSPatrick Sanan char sub_mat_type[256]; 7324b9ad928SBarry Smith 7334b9ad928SBarry Smith PetscFunctionBegin; 734d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject,"Additive Schwarz options"); 7359566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_asm_dm_subdomains","Use DMCreateDomainDecomposition() to define subdomains","PCASMSetDMSubdomains",osm->dm_subdomains,&osm->dm_subdomains,&flg)); 7369566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_asm_blocks","Number of subdomains","PCASMSetTotalSubdomains",osm->n,&blocks,&flg)); 73765db9045SDmitry Karpeev if (flg) { 7389566063dSJacob Faibussowitsch PetscCall(PCASMSetTotalSubdomains(pc,blocks,NULL,NULL)); 739d709ea83SDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 74065db9045SDmitry Karpeev } 7419566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_asm_local_blocks","Number of local subdomains","PCASMSetLocalSubdomains",osm->n_local_true,&blocks,&flg)); 742342c94f9SMatthew G. Knepley if (flg) { 7439566063dSJacob Faibussowitsch PetscCall(PCASMSetLocalSubdomains(pc,blocks,NULL,NULL)); 744342c94f9SMatthew G. Knepley osm->dm_subdomains = PETSC_FALSE; 745342c94f9SMatthew G. Knepley } 7469566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_asm_overlap","Number of grid points overlap","PCASMSetOverlap",osm->overlap,&ovl,&flg)); 74765db9045SDmitry Karpeev if (flg) { 7489566063dSJacob Faibussowitsch PetscCall(PCASMSetOverlap(pc,ovl)); 749d709ea83SDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 75065db9045SDmitry Karpeev } 75190d69ab7SBarry Smith flg = PETSC_FALSE; 7529566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-pc_asm_type","Type of restriction/extension","PCASMSetType",PCASMTypes,(PetscEnum)osm->type,(PetscEnum*)&asmtype,&flg)); 7539566063dSJacob Faibussowitsch if (flg) PetscCall(PCASMSetType(pc,asmtype)); 75412cd4985SMatthew G. Knepley flg = PETSC_FALSE; 7559566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-pc_asm_local_type","Type of local solver composition","PCASMSetLocalType",PCCompositeTypes,(PetscEnum)osm->loctype,(PetscEnum*)&loctype,&flg)); 7569566063dSJacob Faibussowitsch if (flg) PetscCall(PCASMSetLocalType(pc,loctype)); 7579566063dSJacob Faibussowitsch PetscCall(PetscOptionsFList("-pc_asm_sub_mat_type","Subsolve Matrix Type","PCASMSetSubMatType",MatList,NULL,sub_mat_type,256,&flg)); 758*1baa6e33SBarry Smith if (flg) PetscCall(PCASMSetSubMatType(pc,sub_mat_type)); 759d0609cedSBarry Smith PetscOptionsHeadEnd(); 7604b9ad928SBarry Smith PetscFunctionReturn(0); 7614b9ad928SBarry Smith } 7624b9ad928SBarry Smith 7634b9ad928SBarry Smith /*------------------------------------------------------------------------------------*/ 7644b9ad928SBarry Smith 7651e6b0712SBarry Smith static PetscErrorCode PCASMSetLocalSubdomains_ASM(PC pc,PetscInt n,IS is[],IS is_local[]) 7664b9ad928SBarry Smith { 7674b9ad928SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 76892bb6962SLisandro Dalcin PetscInt i; 7694b9ad928SBarry Smith 7704b9ad928SBarry Smith PetscFunctionBegin; 77163a3b9bcSJacob Faibussowitsch PetscCheck(n >= 1,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Each process must have 1 or more blocks, n = %" PetscInt_FMT,n); 7727827d75bSBarry Smith PetscCheck(!pc->setupcalled || (n == osm->n_local_true && !is),PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCASMSetLocalSubdomains() should be called before calling PCSetUp()."); 773e7e72b3dSBarry Smith 7744b9ad928SBarry Smith if (!pc->setupcalled) { 77592bb6962SLisandro Dalcin if (is) { 7769566063dSJacob Faibussowitsch for (i=0; i<n; i++) PetscCall(PetscObjectReference((PetscObject)is[i])); 77792bb6962SLisandro Dalcin } 778832fc9a5SMatthew Knepley if (is_local) { 7799566063dSJacob Faibussowitsch for (i=0; i<n; i++) PetscCall(PetscObjectReference((PetscObject)is_local[i])); 780832fc9a5SMatthew Knepley } 7819566063dSJacob Faibussowitsch PetscCall(PCASMDestroySubdomains(osm->n_local_true,osm->is,osm->is_local)); 7822fa5cd67SKarl Rupp 7834b9ad928SBarry Smith osm->n_local_true = n; 7840a545947SLisandro Dalcin osm->is = NULL; 7850a545947SLisandro Dalcin osm->is_local = NULL; 78692bb6962SLisandro Dalcin if (is) { 7879566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n,&osm->is)); 7882fa5cd67SKarl Rupp for (i=0; i<n; i++) osm->is[i] = is[i]; 7893d03bb48SJed Brown /* Flag indicating that the user has set overlapping subdomains so PCASM should not increase their size. */ 7903d03bb48SJed Brown osm->overlap = -1; 79192bb6962SLisandro Dalcin } 7922b691e39SMatthew Knepley if (is_local) { 7939566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n,&osm->is_local)); 7942fa5cd67SKarl Rupp for (i=0; i<n; i++) osm->is_local[i] = is_local[i]; 795a35b7b57SDmitry Karpeev if (!is) { 7969566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(osm->n_local_true,&osm->is)); 797a35b7b57SDmitry Karpeev for (i=0; i<osm->n_local_true; i++) { 798a35b7b57SDmitry Karpeev if (osm->overlap > 0) { /* With positive overlap, osm->is[i] will be modified */ 7999566063dSJacob Faibussowitsch PetscCall(ISDuplicate(osm->is_local[i],&osm->is[i])); 8009566063dSJacob Faibussowitsch PetscCall(ISCopy(osm->is_local[i],osm->is[i])); 801a35b7b57SDmitry Karpeev } else { 8029566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)osm->is_local[i])); 803a35b7b57SDmitry Karpeev osm->is[i] = osm->is_local[i]; 804a35b7b57SDmitry Karpeev } 805a35b7b57SDmitry Karpeev } 806a35b7b57SDmitry Karpeev } 8072b691e39SMatthew Knepley } 8084b9ad928SBarry Smith } 8094b9ad928SBarry Smith PetscFunctionReturn(0); 8104b9ad928SBarry Smith } 8114b9ad928SBarry Smith 8121e6b0712SBarry Smith static PetscErrorCode PCASMSetTotalSubdomains_ASM(PC pc,PetscInt N,IS *is,IS *is_local) 8134b9ad928SBarry Smith { 8144b9ad928SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 81513f74950SBarry Smith PetscMPIInt rank,size; 81678904715SLisandro Dalcin PetscInt n; 8174b9ad928SBarry Smith 8184b9ad928SBarry Smith PetscFunctionBegin; 81963a3b9bcSJacob Faibussowitsch PetscCheck(N >= 1,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Number of total blocks must be > 0, N = %" PetscInt_FMT,N); 8207827d75bSBarry 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."); 8214b9ad928SBarry Smith 8224b9ad928SBarry Smith /* 823880770e9SJed Brown Split the subdomains equally among all processors 8244b9ad928SBarry Smith */ 8259566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank)); 8269566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size)); 8274b9ad928SBarry Smith n = N/size + ((N % size) > rank); 82863a3b9bcSJacob 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); 8297827d75bSBarry Smith PetscCheck(!pc->setupcalled || n == osm->n_local_true,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"PCASMSetTotalSubdomains() should be called before PCSetUp()."); 8304b9ad928SBarry Smith if (!pc->setupcalled) { 8319566063dSJacob Faibussowitsch PetscCall(PCASMDestroySubdomains(osm->n_local_true,osm->is,osm->is_local)); 8322fa5cd67SKarl Rupp 8334b9ad928SBarry Smith osm->n_local_true = n; 8340a545947SLisandro Dalcin osm->is = NULL; 8350a545947SLisandro Dalcin osm->is_local = NULL; 8364b9ad928SBarry Smith } 8374b9ad928SBarry Smith PetscFunctionReturn(0); 8384b9ad928SBarry Smith } 8394b9ad928SBarry Smith 8401e6b0712SBarry Smith static PetscErrorCode PCASMSetOverlap_ASM(PC pc,PetscInt ovl) 8414b9ad928SBarry Smith { 84292bb6962SLisandro Dalcin PC_ASM *osm = (PC_ASM*)pc->data; 8434b9ad928SBarry Smith 8444b9ad928SBarry Smith PetscFunctionBegin; 8457827d75bSBarry Smith PetscCheck(ovl >= 0,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap value requested"); 8467827d75bSBarry Smith PetscCheck(!pc->setupcalled || ovl == osm->overlap,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCASMSetOverlap() should be called before PCSetUp()."); 8472fa5cd67SKarl Rupp if (!pc->setupcalled) osm->overlap = ovl; 8484b9ad928SBarry Smith PetscFunctionReturn(0); 8494b9ad928SBarry Smith } 8504b9ad928SBarry Smith 8511e6b0712SBarry Smith static PetscErrorCode PCASMSetType_ASM(PC pc,PCASMType type) 8524b9ad928SBarry Smith { 85392bb6962SLisandro Dalcin PC_ASM *osm = (PC_ASM*)pc->data; 8544b9ad928SBarry Smith 8554b9ad928SBarry Smith PetscFunctionBegin; 8564b9ad928SBarry Smith osm->type = type; 857bf108f30SBarry Smith osm->type_set = PETSC_TRUE; 8584b9ad928SBarry Smith PetscFunctionReturn(0); 8594b9ad928SBarry Smith } 8604b9ad928SBarry Smith 861c60c7ad4SBarry Smith static PetscErrorCode PCASMGetType_ASM(PC pc,PCASMType *type) 862c60c7ad4SBarry Smith { 863c60c7ad4SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 864c60c7ad4SBarry Smith 865c60c7ad4SBarry Smith PetscFunctionBegin; 866c60c7ad4SBarry Smith *type = osm->type; 867c60c7ad4SBarry Smith PetscFunctionReturn(0); 868c60c7ad4SBarry Smith } 869c60c7ad4SBarry Smith 87012cd4985SMatthew G. Knepley static PetscErrorCode PCASMSetLocalType_ASM(PC pc, PCCompositeType type) 87112cd4985SMatthew G. Knepley { 87212cd4985SMatthew G. Knepley PC_ASM *osm = (PC_ASM *) pc->data; 87312cd4985SMatthew G. Knepley 87412cd4985SMatthew G. Knepley PetscFunctionBegin; 8757827d75bSBarry 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"); 87612cd4985SMatthew G. Knepley osm->loctype = type; 87712cd4985SMatthew G. Knepley PetscFunctionReturn(0); 87812cd4985SMatthew G. Knepley } 87912cd4985SMatthew G. Knepley 88012cd4985SMatthew G. Knepley static PetscErrorCode PCASMGetLocalType_ASM(PC pc, PCCompositeType *type) 88112cd4985SMatthew G. Knepley { 88212cd4985SMatthew G. Knepley PC_ASM *osm = (PC_ASM *) pc->data; 88312cd4985SMatthew G. Knepley 88412cd4985SMatthew G. Knepley PetscFunctionBegin; 88512cd4985SMatthew G. Knepley *type = osm->loctype; 88612cd4985SMatthew G. Knepley PetscFunctionReturn(0); 88712cd4985SMatthew G. Knepley } 88812cd4985SMatthew G. Knepley 8891e6b0712SBarry Smith static PetscErrorCode PCASMSetSortIndices_ASM(PC pc,PetscBool doSort) 8906ed231c7SMatthew Knepley { 8916ed231c7SMatthew Knepley PC_ASM *osm = (PC_ASM*)pc->data; 8926ed231c7SMatthew Knepley 8936ed231c7SMatthew Knepley PetscFunctionBegin; 8946ed231c7SMatthew Knepley osm->sort_indices = doSort; 8956ed231c7SMatthew Knepley PetscFunctionReturn(0); 8966ed231c7SMatthew Knepley } 8976ed231c7SMatthew Knepley 8981e6b0712SBarry Smith static PetscErrorCode PCASMGetSubKSP_ASM(PC pc,PetscInt *n_local,PetscInt *first_local,KSP **ksp) 8994b9ad928SBarry Smith { 90092bb6962SLisandro Dalcin PC_ASM *osm = (PC_ASM*)pc->data; 9014b9ad928SBarry Smith 9024b9ad928SBarry Smith PetscFunctionBegin; 9037827d75bSBarry 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"); 9044b9ad928SBarry Smith 9052fa5cd67SKarl Rupp if (n_local) *n_local = osm->n_local_true; 90692bb6962SLisandro Dalcin if (first_local) { 9079566063dSJacob Faibussowitsch PetscCallMPI(MPI_Scan(&osm->n_local_true,first_local,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)pc))); 90892bb6962SLisandro Dalcin *first_local -= osm->n_local_true; 90992bb6962SLisandro Dalcin } 910ed155784SPierre Jolivet if (ksp) *ksp = osm->ksp; 9114b9ad928SBarry Smith PetscFunctionReturn(0); 9124b9ad928SBarry Smith } 9134b9ad928SBarry Smith 91480ec0b7dSPatrick Sanan static PetscErrorCode PCASMGetSubMatType_ASM(PC pc,MatType *sub_mat_type) 91580ec0b7dSPatrick Sanan { 91680ec0b7dSPatrick Sanan PC_ASM *osm = (PC_ASM*)pc->data; 91780ec0b7dSPatrick Sanan 91880ec0b7dSPatrick Sanan PetscFunctionBegin; 91980ec0b7dSPatrick Sanan PetscValidHeaderSpecific(pc,PC_CLASSID,1); 92080ec0b7dSPatrick Sanan PetscValidPointer(sub_mat_type,2); 92180ec0b7dSPatrick Sanan *sub_mat_type = osm->sub_mat_type; 92280ec0b7dSPatrick Sanan PetscFunctionReturn(0); 92380ec0b7dSPatrick Sanan } 92480ec0b7dSPatrick Sanan 92580ec0b7dSPatrick Sanan static PetscErrorCode PCASMSetSubMatType_ASM(PC pc,MatType sub_mat_type) 92680ec0b7dSPatrick Sanan { 92780ec0b7dSPatrick Sanan PC_ASM *osm = (PC_ASM*)pc->data; 92880ec0b7dSPatrick Sanan 92980ec0b7dSPatrick Sanan PetscFunctionBegin; 93080ec0b7dSPatrick Sanan PetscValidHeaderSpecific(pc,PC_CLASSID,1); 9319566063dSJacob Faibussowitsch PetscCall(PetscFree(osm->sub_mat_type)); 9329566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(sub_mat_type,(char**)&osm->sub_mat_type)); 93380ec0b7dSPatrick Sanan PetscFunctionReturn(0); 93480ec0b7dSPatrick Sanan } 93580ec0b7dSPatrick Sanan 9364b9ad928SBarry Smith /*@C 9371093a601SBarry Smith PCASMSetLocalSubdomains - Sets the local subdomains (for this processor only) for the additive Schwarz preconditioner. 9384b9ad928SBarry Smith 939d083f849SBarry Smith Collective on pc 9404b9ad928SBarry Smith 9414b9ad928SBarry Smith Input Parameters: 9424b9ad928SBarry Smith + pc - the preconditioner context 9434b9ad928SBarry Smith . n - the number of subdomains for this processor (default value = 1) 9448c03b21aSDmitry Karpeev . is - the index set that defines the subdomains for this processor 9450298fd71SBarry Smith (or NULL for PETSc to determine subdomains) 946f1ee410cSBarry 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 947f1ee410cSBarry Smith (or NULL to not provide these) 9484b9ad928SBarry Smith 949342c94f9SMatthew G. Knepley Options Database Key: 950342c94f9SMatthew G. Knepley To set the total number of subdomain blocks rather than specify the 951342c94f9SMatthew G. Knepley index sets, use the option 952342c94f9SMatthew G. Knepley . -pc_asm_local_blocks <blks> - Sets local blocks 953342c94f9SMatthew G. Knepley 9544b9ad928SBarry Smith Notes: 9551093a601SBarry Smith The IS numbering is in the parallel, global numbering of the vector for both is and is_local 9564b9ad928SBarry Smith 9574b9ad928SBarry Smith By default the ASM preconditioner uses 1 block per processor. 9584b9ad928SBarry Smith 9594b9ad928SBarry Smith Use PCASMSetTotalSubdomains() to set the subdomains for all processors. 9604b9ad928SBarry Smith 961f1ee410cSBarry Smith If is_local is provided and PCASMType is PC_ASM_RESTRICT then the solution only over the is_local region is interpolated 962f1ee410cSBarry Smith back to form the global solution (this is the standard restricted additive Schwarz method) 963f1ee410cSBarry Smith 964f1ee410cSBarry 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 965f1ee410cSBarry Smith no code to handle that case. 966f1ee410cSBarry Smith 9674b9ad928SBarry Smith Level: advanced 9684b9ad928SBarry Smith 969db781477SPatrick Sanan .seealso: `PCASMSetTotalSubdomains()`, `PCASMSetOverlap()`, `PCASMGetSubKSP()`, 970db781477SPatrick Sanan `PCASMCreateSubdomains2D()`, `PCASMGetLocalSubdomains()`, `PCASMType`, `PCASMSetType()` 9714b9ad928SBarry Smith @*/ 9727087cfbeSBarry Smith PetscErrorCode PCASMSetLocalSubdomains(PC pc,PetscInt n,IS is[],IS is_local[]) 9734b9ad928SBarry Smith { 9744b9ad928SBarry Smith PetscFunctionBegin; 9750700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 976cac4c232SBarry Smith PetscTryMethod(pc,"PCASMSetLocalSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,n,is,is_local)); 9774b9ad928SBarry Smith PetscFunctionReturn(0); 9784b9ad928SBarry Smith } 9794b9ad928SBarry Smith 9804b9ad928SBarry Smith /*@C 981feb221f8SDmitry Karpeev PCASMSetTotalSubdomains - Sets the subdomains for all processors for the 9824b9ad928SBarry Smith additive Schwarz preconditioner. Either all or no processors in the 9834b9ad928SBarry Smith PC communicator must call this routine, with the same index sets. 9844b9ad928SBarry Smith 985d083f849SBarry Smith Collective on pc 9864b9ad928SBarry Smith 9874b9ad928SBarry Smith Input Parameters: 9884b9ad928SBarry Smith + pc - the preconditioner context 989feb221f8SDmitry Karpeev . N - the number of subdomains for all processors 990feb221f8SDmitry Karpeev . is - the index sets that define the subdomains for all processors 991dfaa0c5fSBarry Smith (or NULL to ask PETSc to determine the subdomains) 9922b691e39SMatthew Knepley - is_local - the index sets that define the local part of the subdomains for this processor 993f1ee410cSBarry Smith (or NULL to not provide this information) 9944b9ad928SBarry Smith 9954b9ad928SBarry Smith Options Database Key: 9964b9ad928SBarry Smith To set the total number of subdomain blocks rather than specify the 9974b9ad928SBarry Smith index sets, use the option 9984b9ad928SBarry Smith . -pc_asm_blocks <blks> - Sets total blocks 9994b9ad928SBarry Smith 10004b9ad928SBarry Smith Notes: 1001f1ee410cSBarry Smith Currently you cannot use this to set the actual subdomains with the argument is or is_local. 10024b9ad928SBarry Smith 10034b9ad928SBarry Smith By default the ASM preconditioner uses 1 block per processor. 10044b9ad928SBarry Smith 10054b9ad928SBarry Smith These index sets cannot be destroyed until after completion of the 10064b9ad928SBarry Smith linear solves for which the ASM preconditioner is being used. 10074b9ad928SBarry Smith 10084b9ad928SBarry Smith Use PCASMSetLocalSubdomains() to set local subdomains. 10094b9ad928SBarry Smith 10101093a601SBarry Smith The IS numbering is in the parallel, global numbering of the vector for both is and is_local 10111093a601SBarry Smith 10124b9ad928SBarry Smith Level: advanced 10134b9ad928SBarry Smith 1014db781477SPatrick Sanan .seealso: `PCASMSetLocalSubdomains()`, `PCASMSetOverlap()`, `PCASMGetSubKSP()`, 1015db781477SPatrick Sanan `PCASMCreateSubdomains2D()` 10164b9ad928SBarry Smith @*/ 10177087cfbeSBarry Smith PetscErrorCode PCASMSetTotalSubdomains(PC pc,PetscInt N,IS is[],IS is_local[]) 10184b9ad928SBarry Smith { 10194b9ad928SBarry Smith PetscFunctionBegin; 10200700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1021cac4c232SBarry Smith PetscTryMethod(pc,"PCASMSetTotalSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,N,is,is_local)); 10224b9ad928SBarry Smith PetscFunctionReturn(0); 10234b9ad928SBarry Smith } 10244b9ad928SBarry Smith 10254b9ad928SBarry Smith /*@ 10264b9ad928SBarry Smith PCASMSetOverlap - Sets the overlap between a pair of subdomains for the 10274b9ad928SBarry Smith additive Schwarz preconditioner. Either all or no processors in the 1028f1ee410cSBarry Smith PC communicator must call this routine. 10294b9ad928SBarry Smith 1030d083f849SBarry Smith Logically Collective on pc 10314b9ad928SBarry Smith 10324b9ad928SBarry Smith Input Parameters: 10334b9ad928SBarry Smith + pc - the preconditioner context 10344b9ad928SBarry Smith - ovl - the amount of overlap between subdomains (ovl >= 0, default value = 1) 10354b9ad928SBarry Smith 10364b9ad928SBarry Smith Options Database Key: 10374b9ad928SBarry Smith . -pc_asm_overlap <ovl> - Sets overlap 10384b9ad928SBarry Smith 10394b9ad928SBarry Smith Notes: 10404b9ad928SBarry Smith By default the ASM preconditioner uses 1 block per processor. To use 10414b9ad928SBarry Smith multiple blocks per perocessor, see PCASMSetTotalSubdomains() and 10424b9ad928SBarry Smith PCASMSetLocalSubdomains() (and the option -pc_asm_blocks <blks>). 10434b9ad928SBarry Smith 10444b9ad928SBarry Smith The overlap defaults to 1, so if one desires that no additional 10454b9ad928SBarry Smith overlap be computed beyond what may have been set with a call to 10464b9ad928SBarry Smith PCASMSetTotalSubdomains() or PCASMSetLocalSubdomains(), then ovl 10474b9ad928SBarry Smith must be set to be 0. In particular, if one does not explicitly set 10484b9ad928SBarry Smith the subdomains an application code, then all overlap would be computed 10494b9ad928SBarry Smith internally by PETSc, and using an overlap of 0 would result in an ASM 10504b9ad928SBarry Smith variant that is equivalent to the block Jacobi preconditioner. 10514b9ad928SBarry Smith 1052f1ee410cSBarry Smith The default algorithm used by PETSc to increase overlap is fast, but not scalable, 1053f1ee410cSBarry Smith use the option -mat_increase_overlap_scalable when the problem and number of processes is large. 1054f1ee410cSBarry Smith 10554b9ad928SBarry Smith Note that one can define initial index sets with any overlap via 1056f1ee410cSBarry Smith PCASMSetLocalSubdomains(); the routine 10574b9ad928SBarry Smith PCASMSetOverlap() merely allows PETSc to extend that overlap further 10584b9ad928SBarry Smith if desired. 10594b9ad928SBarry Smith 10604b9ad928SBarry Smith Level: intermediate 10614b9ad928SBarry Smith 1062db781477SPatrick Sanan .seealso: `PCASMSetTotalSubdomains()`, `PCASMSetLocalSubdomains()`, `PCASMGetSubKSP()`, 1063db781477SPatrick Sanan `PCASMCreateSubdomains2D()`, `PCASMGetLocalSubdomains()`, `MatIncreaseOverlap()` 10644b9ad928SBarry Smith @*/ 10657087cfbeSBarry Smith PetscErrorCode PCASMSetOverlap(PC pc,PetscInt ovl) 10664b9ad928SBarry Smith { 10674b9ad928SBarry Smith PetscFunctionBegin; 10680700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1069c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,ovl,2); 1070cac4c232SBarry Smith PetscTryMethod(pc,"PCASMSetOverlap_C",(PC,PetscInt),(pc,ovl)); 10714b9ad928SBarry Smith PetscFunctionReturn(0); 10724b9ad928SBarry Smith } 10734b9ad928SBarry Smith 10744b9ad928SBarry Smith /*@ 10754b9ad928SBarry Smith PCASMSetType - Sets the type of restriction and interpolation used 10764b9ad928SBarry Smith for local problems in the additive Schwarz method. 10774b9ad928SBarry Smith 1078d083f849SBarry Smith Logically Collective on pc 10794b9ad928SBarry Smith 10804b9ad928SBarry Smith Input Parameters: 10814b9ad928SBarry Smith + pc - the preconditioner context 10824b9ad928SBarry Smith - type - variant of ASM, one of 10834b9ad928SBarry Smith .vb 10844b9ad928SBarry Smith PC_ASM_BASIC - full interpolation and restriction 108582b5ce2aSStefano Zampini PC_ASM_RESTRICT - full restriction, local processor interpolation (default) 10864b9ad928SBarry Smith PC_ASM_INTERPOLATE - full interpolation, local processor restriction 10874b9ad928SBarry Smith PC_ASM_NONE - local processor restriction and interpolation 10884b9ad928SBarry Smith .ve 10894b9ad928SBarry Smith 10904b9ad928SBarry Smith Options Database Key: 10914b9ad928SBarry Smith . -pc_asm_type [basic,restrict,interpolate,none] - Sets ASM type 10924b9ad928SBarry Smith 109395452b02SPatrick Sanan Notes: 109495452b02SPatrick Sanan if the is_local arguments are passed to PCASMSetLocalSubdomains() then they are used when PC_ASM_RESTRICT has been selected 1095f1ee410cSBarry Smith to limit the local processor interpolation 1096f1ee410cSBarry Smith 10974b9ad928SBarry Smith Level: intermediate 10984b9ad928SBarry Smith 1099db781477SPatrick Sanan .seealso: `PCASMSetTotalSubdomains()`, `PCASMSetTotalSubdomains()`, `PCASMGetSubKSP()`, 1100db781477SPatrick Sanan `PCASMCreateSubdomains2D()`, `PCASMType`, `PCASMSetLocalType()`, `PCASMGetLocalType()` 11014b9ad928SBarry Smith @*/ 11027087cfbeSBarry Smith PetscErrorCode PCASMSetType(PC pc,PCASMType type) 11034b9ad928SBarry Smith { 11044b9ad928SBarry Smith PetscFunctionBegin; 11050700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1106c5eb9154SBarry Smith PetscValidLogicalCollectiveEnum(pc,type,2); 1107cac4c232SBarry Smith PetscTryMethod(pc,"PCASMSetType_C",(PC,PCASMType),(pc,type)); 11084b9ad928SBarry Smith PetscFunctionReturn(0); 11094b9ad928SBarry Smith } 11104b9ad928SBarry Smith 1111c60c7ad4SBarry Smith /*@ 1112c60c7ad4SBarry Smith PCASMGetType - Gets the type of restriction and interpolation used 1113c60c7ad4SBarry Smith for local problems in the additive Schwarz method. 1114c60c7ad4SBarry Smith 1115d083f849SBarry Smith Logically Collective on pc 1116c60c7ad4SBarry Smith 1117c60c7ad4SBarry Smith Input Parameter: 1118c60c7ad4SBarry Smith . pc - the preconditioner context 1119c60c7ad4SBarry Smith 1120c60c7ad4SBarry Smith Output Parameter: 1121c60c7ad4SBarry Smith . type - variant of ASM, one of 1122c60c7ad4SBarry Smith 1123c60c7ad4SBarry Smith .vb 1124c60c7ad4SBarry Smith PC_ASM_BASIC - full interpolation and restriction 1125c60c7ad4SBarry Smith PC_ASM_RESTRICT - full restriction, local processor interpolation 1126c60c7ad4SBarry Smith PC_ASM_INTERPOLATE - full interpolation, local processor restriction 1127c60c7ad4SBarry Smith PC_ASM_NONE - local processor restriction and interpolation 1128c60c7ad4SBarry Smith .ve 1129c60c7ad4SBarry Smith 1130c60c7ad4SBarry Smith Options Database Key: 1131c60c7ad4SBarry Smith . -pc_asm_type [basic,restrict,interpolate,none] - Sets ASM type 1132c60c7ad4SBarry Smith 1133c60c7ad4SBarry Smith Level: intermediate 1134c60c7ad4SBarry Smith 1135db781477SPatrick Sanan .seealso: `PCASMSetTotalSubdomains()`, `PCASMSetTotalSubdomains()`, `PCASMGetSubKSP()`, 1136db781477SPatrick Sanan `PCASMCreateSubdomains2D()`, `PCASMType`, `PCASMSetType()`, `PCASMSetLocalType()`, `PCASMGetLocalType()` 1137c60c7ad4SBarry Smith @*/ 1138c60c7ad4SBarry Smith PetscErrorCode PCASMGetType(PC pc,PCASMType *type) 1139c60c7ad4SBarry Smith { 1140c60c7ad4SBarry Smith PetscFunctionBegin; 1141c60c7ad4SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1142cac4c232SBarry Smith PetscUseMethod(pc,"PCASMGetType_C",(PC,PCASMType*),(pc,type)); 1143c60c7ad4SBarry Smith PetscFunctionReturn(0); 1144c60c7ad4SBarry Smith } 1145c60c7ad4SBarry Smith 114612cd4985SMatthew G. Knepley /*@ 114712cd4985SMatthew G. Knepley PCASMSetLocalType - Sets the type of composition used for local problems in the additive Schwarz method. 114812cd4985SMatthew G. Knepley 1149d083f849SBarry Smith Logically Collective on pc 115012cd4985SMatthew G. Knepley 115112cd4985SMatthew G. Knepley Input Parameters: 115212cd4985SMatthew G. Knepley + pc - the preconditioner context 115312cd4985SMatthew G. Knepley - type - type of composition, one of 115412cd4985SMatthew G. Knepley .vb 115512cd4985SMatthew G. Knepley PC_COMPOSITE_ADDITIVE - local additive combination 115612cd4985SMatthew G. Knepley PC_COMPOSITE_MULTIPLICATIVE - local multiplicative combination 115712cd4985SMatthew G. Knepley .ve 115812cd4985SMatthew G. Knepley 115912cd4985SMatthew G. Knepley Options Database Key: 116012cd4985SMatthew G. Knepley . -pc_asm_local_type [additive,multiplicative] - Sets local solver composition type 116112cd4985SMatthew G. Knepley 116212cd4985SMatthew G. Knepley Level: intermediate 116312cd4985SMatthew G. Knepley 1164db781477SPatrick Sanan .seealso: `PCASMSetType()`, `PCASMGetType()`, `PCASMGetLocalType()`, `PCASM`, `PCASMType`, `PCASMSetType()`, `PCASMGetType()`, `PCCompositeType` 116512cd4985SMatthew G. Knepley @*/ 116612cd4985SMatthew G. Knepley PetscErrorCode PCASMSetLocalType(PC pc, PCCompositeType type) 116712cd4985SMatthew G. Knepley { 116812cd4985SMatthew G. Knepley PetscFunctionBegin; 116912cd4985SMatthew G. Knepley PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 117012cd4985SMatthew G. Knepley PetscValidLogicalCollectiveEnum(pc, type, 2); 1171cac4c232SBarry Smith PetscTryMethod(pc, "PCASMSetLocalType_C", (PC, PCCompositeType), (pc, type)); 117212cd4985SMatthew G. Knepley PetscFunctionReturn(0); 117312cd4985SMatthew G. Knepley } 117412cd4985SMatthew G. Knepley 117512cd4985SMatthew G. Knepley /*@ 117612cd4985SMatthew G. Knepley PCASMGetLocalType - Gets the type of composition used for local problems in the additive Schwarz method. 117712cd4985SMatthew G. Knepley 1178d083f849SBarry Smith Logically Collective on pc 117912cd4985SMatthew G. Knepley 118012cd4985SMatthew G. Knepley Input Parameter: 118112cd4985SMatthew G. Knepley . pc - the preconditioner context 118212cd4985SMatthew G. Knepley 118312cd4985SMatthew G. Knepley Output Parameter: 118412cd4985SMatthew G. Knepley . type - type of composition, one of 118512cd4985SMatthew G. Knepley .vb 118612cd4985SMatthew G. Knepley PC_COMPOSITE_ADDITIVE - local additive combination 118712cd4985SMatthew G. Knepley PC_COMPOSITE_MULTIPLICATIVE - local multiplicative combination 118812cd4985SMatthew G. Knepley .ve 118912cd4985SMatthew G. Knepley 119012cd4985SMatthew G. Knepley Options Database Key: 119112cd4985SMatthew G. Knepley . -pc_asm_local_type [additive,multiplicative] - Sets local solver composition type 119212cd4985SMatthew G. Knepley 119312cd4985SMatthew G. Knepley Level: intermediate 119412cd4985SMatthew G. Knepley 1195db781477SPatrick Sanan .seealso: `PCASMSetType()`, `PCASMGetType()`, `PCASMSetLocalType()`, `PCASMCreate()`, `PCASMType`, `PCASMSetType()`, `PCASMGetType()`, `PCCompositeType` 119612cd4985SMatthew G. Knepley @*/ 119712cd4985SMatthew G. Knepley PetscErrorCode PCASMGetLocalType(PC pc, PCCompositeType *type) 119812cd4985SMatthew G. Knepley { 119912cd4985SMatthew G. Knepley PetscFunctionBegin; 120012cd4985SMatthew G. Knepley PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 120112cd4985SMatthew G. Knepley PetscValidPointer(type, 2); 1202cac4c232SBarry Smith PetscUseMethod(pc, "PCASMGetLocalType_C", (PC, PCCompositeType *), (pc, type)); 120312cd4985SMatthew G. Knepley PetscFunctionReturn(0); 120412cd4985SMatthew G. Knepley } 120512cd4985SMatthew G. Knepley 12066ed231c7SMatthew Knepley /*@ 12076ed231c7SMatthew Knepley PCASMSetSortIndices - Determines whether subdomain indices are sorted. 12086ed231c7SMatthew Knepley 1209d083f849SBarry Smith Logically Collective on pc 12106ed231c7SMatthew Knepley 12116ed231c7SMatthew Knepley Input Parameters: 12126ed231c7SMatthew Knepley + pc - the preconditioner context 12136ed231c7SMatthew Knepley - doSort - sort the subdomain indices 12146ed231c7SMatthew Knepley 12156ed231c7SMatthew Knepley Level: intermediate 12166ed231c7SMatthew Knepley 1217db781477SPatrick Sanan .seealso: `PCASMSetLocalSubdomains()`, `PCASMSetTotalSubdomains()`, `PCASMGetSubKSP()`, 1218db781477SPatrick Sanan `PCASMCreateSubdomains2D()` 12196ed231c7SMatthew Knepley @*/ 12207087cfbeSBarry Smith PetscErrorCode PCASMSetSortIndices(PC pc,PetscBool doSort) 12216ed231c7SMatthew Knepley { 12226ed231c7SMatthew Knepley PetscFunctionBegin; 12230700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1224acfcf0e5SJed Brown PetscValidLogicalCollectiveBool(pc,doSort,2); 1225cac4c232SBarry Smith PetscTryMethod(pc,"PCASMSetSortIndices_C",(PC,PetscBool),(pc,doSort)); 12266ed231c7SMatthew Knepley PetscFunctionReturn(0); 12276ed231c7SMatthew Knepley } 12286ed231c7SMatthew Knepley 12294b9ad928SBarry Smith /*@C 12304b9ad928SBarry Smith PCASMGetSubKSP - Gets the local KSP contexts for all blocks on 12314b9ad928SBarry Smith this processor. 12324b9ad928SBarry Smith 1233d083f849SBarry Smith Collective on pc iff first_local is requested 12344b9ad928SBarry Smith 12354b9ad928SBarry Smith Input Parameter: 12364b9ad928SBarry Smith . pc - the preconditioner context 12374b9ad928SBarry Smith 12384b9ad928SBarry Smith Output Parameters: 12390298fd71SBarry Smith + n_local - the number of blocks on this processor or NULL 12400298fd71SBarry Smith . first_local - the global number of the first block on this processor or NULL, 12410298fd71SBarry Smith all processors must request or all must pass NULL 12424b9ad928SBarry Smith - ksp - the array of KSP contexts 12434b9ad928SBarry Smith 12444b9ad928SBarry Smith Note: 1245d29017ddSJed Brown After PCASMGetSubKSP() the array of KSPes is not to be freed. 12464b9ad928SBarry Smith 12474b9ad928SBarry Smith You must call KSPSetUp() before calling PCASMGetSubKSP(). 12484b9ad928SBarry Smith 1249d29017ddSJed Brown Fortran note: 12502bf68e3eSBarry 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. 1251d29017ddSJed Brown 12524b9ad928SBarry Smith Level: advanced 12534b9ad928SBarry Smith 1254db781477SPatrick Sanan .seealso: `PCASMSetTotalSubdomains()`, `PCASMSetTotalSubdomains()`, `PCASMSetOverlap()`, 1255db781477SPatrick Sanan `PCASMCreateSubdomains2D()`, 12564b9ad928SBarry Smith @*/ 12577087cfbeSBarry Smith PetscErrorCode PCASMGetSubKSP(PC pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp[]) 12584b9ad928SBarry Smith { 12594b9ad928SBarry Smith PetscFunctionBegin; 12600700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1261cac4c232SBarry Smith PetscUseMethod(pc,"PCASMGetSubKSP_C",(PC,PetscInt*,PetscInt*,KSP **),(pc,n_local,first_local,ksp)); 12624b9ad928SBarry Smith PetscFunctionReturn(0); 12634b9ad928SBarry Smith } 12644b9ad928SBarry Smith 12654b9ad928SBarry Smith /* -------------------------------------------------------------------------------------*/ 12664b9ad928SBarry Smith /*MC 12673b09bd56SBarry Smith PCASM - Use the (restricted) additive Schwarz method, each block is (approximately) solved with 12684b9ad928SBarry Smith its own KSP object. 12694b9ad928SBarry Smith 12704b9ad928SBarry Smith Options Database Keys: 127149517cdeSBarry Smith + -pc_asm_blocks <blks> - Sets total blocks 12724b9ad928SBarry Smith . -pc_asm_overlap <ovl> - Sets overlap 1273f1ee410cSBarry Smith . -pc_asm_type [basic,restrict,interpolate,none] - Sets ASM type, default is restrict 1274f1ee410cSBarry Smith - -pc_asm_local_type [additive, multiplicative] - Sets ASM type, default is additive 12754b9ad928SBarry Smith 12763b09bd56SBarry Smith IMPORTANT: If you run with, for example, 3 blocks on 1 processor or 3 blocks on 3 processors you 12773b09bd56SBarry Smith will get a different convergence rate due to the default option of -pc_asm_type restrict. Use 12783b09bd56SBarry Smith -pc_asm_type basic to use the standard ASM. 12793b09bd56SBarry Smith 128095452b02SPatrick Sanan Notes: 128195452b02SPatrick Sanan Each processor can have one or more blocks, but a block cannot be shared by more 1282f1ee410cSBarry Smith than one processor. Use PCGASM for subdomains shared by multiple processes. Defaults to one block per processor. 12834b9ad928SBarry Smith 12843b09bd56SBarry Smith To set options on the solvers for each block append -sub_ to all the KSP, and PC 1285d7ee0231SBarry Smith options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly 12864b9ad928SBarry Smith 1287a8c7a070SBarry Smith To set the options on the solvers separate for each block call PCASMGetSubKSP() 12884b9ad928SBarry Smith and set the options directly on the resulting KSP object (you can access its PC 12894b9ad928SBarry Smith with KSPGetPC()) 12904b9ad928SBarry Smith 12914b9ad928SBarry Smith Level: beginner 12924b9ad928SBarry Smith 1293c582cd25SBarry Smith References: 1294606c0280SSatish Balay + * - M Dryja, OB Widlund, An additive variant of the Schwarz alternating method for the case of many subregions 129596a0c994SBarry Smith Courant Institute, New York University Technical report 1296606c0280SSatish Balay - * - Barry Smith, Petter Bjorstad, and William Gropp, Domain Decompositions: Parallel Multilevel Methods for Elliptic Partial Differential Equations, 129796a0c994SBarry Smith Cambridge University Press. 1298c582cd25SBarry Smith 1299db781477SPatrick Sanan .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, 1300db781477SPatrick Sanan `PCBJACOBI`, `PCASMGetSubKSP()`, `PCASMSetLocalSubdomains()`, `PCASMType`, `PCASMGetType()`, `PCASMSetLocalType()`, `PCASMGetLocalType()` 1301db781477SPatrick Sanan `PCASMSetTotalSubdomains()`, `PCSetModifySubMatrices()`, `PCASMSetOverlap()`, `PCASMSetType()`, `PCCompositeType` 1302e09e08ccSBarry Smith 13034b9ad928SBarry Smith M*/ 13044b9ad928SBarry Smith 13058cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_ASM(PC pc) 13064b9ad928SBarry Smith { 13074b9ad928SBarry Smith PC_ASM *osm; 13084b9ad928SBarry Smith 13094b9ad928SBarry Smith PetscFunctionBegin; 13109566063dSJacob Faibussowitsch PetscCall(PetscNewLog(pc,&osm)); 13112fa5cd67SKarl Rupp 13124b9ad928SBarry Smith osm->n = PETSC_DECIDE; 13134b9ad928SBarry Smith osm->n_local = 0; 13142b837212SDmitry Karpeev osm->n_local_true = PETSC_DECIDE; 13154b9ad928SBarry Smith osm->overlap = 1; 13160a545947SLisandro Dalcin osm->ksp = NULL; 13170a545947SLisandro Dalcin osm->restriction = NULL; 13180a545947SLisandro Dalcin osm->lprolongation = NULL; 13190a545947SLisandro Dalcin osm->lrestriction = NULL; 13200a545947SLisandro Dalcin osm->x = NULL; 13210a545947SLisandro Dalcin osm->y = NULL; 13220a545947SLisandro Dalcin osm->is = NULL; 13230a545947SLisandro Dalcin osm->is_local = NULL; 13240a545947SLisandro Dalcin osm->mat = NULL; 13250a545947SLisandro Dalcin osm->pmat = NULL; 13264b9ad928SBarry Smith osm->type = PC_ASM_RESTRICT; 132712cd4985SMatthew G. Knepley osm->loctype = PC_COMPOSITE_ADDITIVE; 13286ed231c7SMatthew Knepley osm->sort_indices = PETSC_TRUE; 1329d709ea83SDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 133080ec0b7dSPatrick Sanan osm->sub_mat_type = NULL; 13314b9ad928SBarry Smith 133292bb6962SLisandro Dalcin pc->data = (void*)osm; 13334b9ad928SBarry Smith pc->ops->apply = PCApply_ASM; 133448e38a8aSPierre Jolivet pc->ops->matapply = PCMatApply_ASM; 13354b9ad928SBarry Smith pc->ops->applytranspose = PCApplyTranspose_ASM; 13364b9ad928SBarry Smith pc->ops->setup = PCSetUp_ASM; 1337e91c6855SBarry Smith pc->ops->reset = PCReset_ASM; 13384b9ad928SBarry Smith pc->ops->destroy = PCDestroy_ASM; 13394b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_ASM; 13404b9ad928SBarry Smith pc->ops->setuponblocks = PCSetUpOnBlocks_ASM; 13414b9ad928SBarry Smith pc->ops->view = PCView_ASM; 13420a545947SLisandro Dalcin pc->ops->applyrichardson = NULL; 13434b9ad928SBarry Smith 13449566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCASMSetLocalSubdomains_C",PCASMSetLocalSubdomains_ASM)); 13459566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCASMSetTotalSubdomains_C",PCASMSetTotalSubdomains_ASM)); 13469566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCASMSetOverlap_C",PCASMSetOverlap_ASM)); 13479566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCASMSetType_C",PCASMSetType_ASM)); 13489566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCASMGetType_C",PCASMGetType_ASM)); 13499566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCASMSetLocalType_C",PCASMSetLocalType_ASM)); 13509566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCASMGetLocalType_C",PCASMGetLocalType_ASM)); 13519566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCASMSetSortIndices_C",PCASMSetSortIndices_ASM)); 13529566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCASMGetSubKSP_C",PCASMGetSubKSP_ASM)); 13539566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCASMGetSubMatType_C",PCASMGetSubMatType_ASM)); 13549566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCASMSetSubMatType_C",PCASMSetSubMatType_ASM)); 13554b9ad928SBarry Smith PetscFunctionReturn(0); 13564b9ad928SBarry Smith } 13574b9ad928SBarry Smith 135892bb6962SLisandro Dalcin /*@C 135992bb6962SLisandro Dalcin PCASMCreateSubdomains - Creates the index sets for the overlapping Schwarz 13603b8d980fSPierre Jolivet preconditioner for any problem on a general grid. 136192bb6962SLisandro Dalcin 136292bb6962SLisandro Dalcin Collective 136392bb6962SLisandro Dalcin 136492bb6962SLisandro Dalcin Input Parameters: 136592bb6962SLisandro Dalcin + A - The global matrix operator 136692bb6962SLisandro Dalcin - n - the number of local blocks 136792bb6962SLisandro Dalcin 136892bb6962SLisandro Dalcin Output Parameters: 136992bb6962SLisandro Dalcin . outis - the array of index sets defining the subdomains 137092bb6962SLisandro Dalcin 137192bb6962SLisandro Dalcin Level: advanced 137292bb6962SLisandro Dalcin 13737d6bfa3bSBarry Smith Note: this generates nonoverlapping subdomains; the PCASM will generate the overlap 13747d6bfa3bSBarry Smith from these if you use PCASMSetLocalSubdomains() 13757d6bfa3bSBarry Smith 13767d6bfa3bSBarry Smith In the Fortran version you must provide the array outis[] already allocated of length n. 13777d6bfa3bSBarry Smith 1378db781477SPatrick Sanan .seealso: `PCASMSetLocalSubdomains()`, `PCASMDestroySubdomains()` 137992bb6962SLisandro Dalcin @*/ 13807087cfbeSBarry Smith PetscErrorCode PCASMCreateSubdomains(Mat A, PetscInt n, IS* outis[]) 138192bb6962SLisandro Dalcin { 138292bb6962SLisandro Dalcin MatPartitioning mpart; 138392bb6962SLisandro Dalcin const char *prefix; 138492bb6962SLisandro Dalcin PetscInt i,j,rstart,rend,bs; 1385976e8c5aSLisandro Dalcin PetscBool hasop, isbaij = PETSC_FALSE,foundpart = PETSC_FALSE; 13860298fd71SBarry Smith Mat Ad = NULL, adj; 138792bb6962SLisandro Dalcin IS ispart,isnumb,*is; 138892bb6962SLisandro Dalcin 138992bb6962SLisandro Dalcin PetscFunctionBegin; 13900700a824SBarry Smith PetscValidHeaderSpecific(A,MAT_CLASSID,1); 139192bb6962SLisandro Dalcin PetscValidPointer(outis,3); 139263a3b9bcSJacob Faibussowitsch PetscCheck(n >= 1,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of local blocks must be > 0, n = %" PetscInt_FMT,n); 139392bb6962SLisandro Dalcin 139492bb6962SLisandro Dalcin /* Get prefix, row distribution, and block size */ 13959566063dSJacob Faibussowitsch PetscCall(MatGetOptionsPrefix(A,&prefix)); 13969566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A,&rstart,&rend)); 13979566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(A,&bs)); 139863a3b9bcSJacob 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); 139965e19b50SBarry Smith 140092bb6962SLisandro Dalcin /* Get diagonal block from matrix if possible */ 14019566063dSJacob Faibussowitsch PetscCall(MatHasOperation(A,MATOP_GET_DIAGONAL_BLOCK,&hasop)); 1402976e8c5aSLisandro Dalcin if (hasop) { 14039566063dSJacob Faibussowitsch PetscCall(MatGetDiagonalBlock(A,&Ad)); 140492bb6962SLisandro Dalcin } 140592bb6962SLisandro Dalcin if (Ad) { 14069566063dSJacob Faibussowitsch PetscCall(PetscObjectBaseTypeCompare((PetscObject)Ad,MATSEQBAIJ,&isbaij)); 14079566063dSJacob Faibussowitsch if (!isbaij) PetscCall(PetscObjectBaseTypeCompare((PetscObject)Ad,MATSEQSBAIJ,&isbaij)); 140892bb6962SLisandro Dalcin } 140992bb6962SLisandro Dalcin if (Ad && n > 1) { 1410ace3abfcSBarry Smith PetscBool match,done; 141192bb6962SLisandro Dalcin /* Try to setup a good matrix partitioning if available */ 14129566063dSJacob Faibussowitsch PetscCall(MatPartitioningCreate(PETSC_COMM_SELF,&mpart)); 14139566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix)); 14149566063dSJacob Faibussowitsch PetscCall(MatPartitioningSetFromOptions(mpart)); 14159566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGCURRENT,&match)); 141692bb6962SLisandro Dalcin if (!match) { 14179566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGSQUARE,&match)); 141892bb6962SLisandro Dalcin } 141992bb6962SLisandro Dalcin if (!match) { /* assume a "good" partitioner is available */ 14201a83f524SJed Brown PetscInt na; 14211a83f524SJed Brown const PetscInt *ia,*ja; 14229566063dSJacob Faibussowitsch PetscCall(MatGetRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done)); 142392bb6962SLisandro Dalcin if (done) { 142492bb6962SLisandro Dalcin /* Build adjacency matrix by hand. Unfortunately a call to 142592bb6962SLisandro Dalcin MatConvert(Ad,MATMPIADJ,MAT_INITIAL_MATRIX,&adj) will 142692bb6962SLisandro Dalcin remove the block-aij structure and we cannot expect 142792bb6962SLisandro Dalcin MatPartitioning to split vertices as we need */ 14280a545947SLisandro Dalcin PetscInt i,j,len,nnz,cnt,*iia=NULL,*jja=NULL; 14291a83f524SJed Brown const PetscInt *row; 143092bb6962SLisandro Dalcin nnz = 0; 143192bb6962SLisandro Dalcin for (i=0; i<na; i++) { /* count number of nonzeros */ 143292bb6962SLisandro Dalcin len = ia[i+1] - ia[i]; 143392bb6962SLisandro Dalcin row = ja + ia[i]; 143492bb6962SLisandro Dalcin for (j=0; j<len; j++) { 143592bb6962SLisandro Dalcin if (row[j] == i) { /* don't count diagonal */ 143692bb6962SLisandro Dalcin len--; break; 143792bb6962SLisandro Dalcin } 143892bb6962SLisandro Dalcin } 143992bb6962SLisandro Dalcin nnz += len; 144092bb6962SLisandro Dalcin } 14419566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(na+1,&iia)); 14429566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nnz,&jja)); 144392bb6962SLisandro Dalcin nnz = 0; 144492bb6962SLisandro Dalcin iia[0] = 0; 144592bb6962SLisandro Dalcin for (i=0; i<na; i++) { /* fill adjacency */ 144692bb6962SLisandro Dalcin cnt = 0; 144792bb6962SLisandro Dalcin len = ia[i+1] - ia[i]; 144892bb6962SLisandro Dalcin row = ja + ia[i]; 144992bb6962SLisandro Dalcin for (j=0; j<len; j++) { 145092bb6962SLisandro Dalcin if (row[j] != i) { /* if not diagonal */ 145192bb6962SLisandro Dalcin jja[nnz+cnt++] = row[j]; 145292bb6962SLisandro Dalcin } 145392bb6962SLisandro Dalcin } 145492bb6962SLisandro Dalcin nnz += cnt; 145592bb6962SLisandro Dalcin iia[i+1] = nnz; 145692bb6962SLisandro Dalcin } 145792bb6962SLisandro Dalcin /* Partitioning of the adjacency matrix */ 14589566063dSJacob Faibussowitsch PetscCall(MatCreateMPIAdj(PETSC_COMM_SELF,na,na,iia,jja,NULL,&adj)); 14599566063dSJacob Faibussowitsch PetscCall(MatPartitioningSetAdjacency(mpart,adj)); 14609566063dSJacob Faibussowitsch PetscCall(MatPartitioningSetNParts(mpart,n)); 14619566063dSJacob Faibussowitsch PetscCall(MatPartitioningApply(mpart,&ispart)); 14629566063dSJacob Faibussowitsch PetscCall(ISPartitioningToNumbering(ispart,&isnumb)); 14639566063dSJacob Faibussowitsch PetscCall(MatDestroy(&adj)); 146492bb6962SLisandro Dalcin foundpart = PETSC_TRUE; 146592bb6962SLisandro Dalcin } 14669566063dSJacob Faibussowitsch PetscCall(MatRestoreRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done)); 146792bb6962SLisandro Dalcin } 14689566063dSJacob Faibussowitsch PetscCall(MatPartitioningDestroy(&mpart)); 146992bb6962SLisandro Dalcin } 147092bb6962SLisandro Dalcin 14719566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n,&is)); 147292bb6962SLisandro Dalcin *outis = is; 147392bb6962SLisandro Dalcin 147492bb6962SLisandro Dalcin if (!foundpart) { 147592bb6962SLisandro Dalcin 147692bb6962SLisandro Dalcin /* Partitioning by contiguous chunks of rows */ 147792bb6962SLisandro Dalcin 147892bb6962SLisandro Dalcin PetscInt mbs = (rend-rstart)/bs; 147992bb6962SLisandro Dalcin PetscInt start = rstart; 148092bb6962SLisandro Dalcin for (i=0; i<n; i++) { 148192bb6962SLisandro Dalcin PetscInt count = (mbs/n + ((mbs % n) > i)) * bs; 14829566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF,count,start,1,&is[i])); 148392bb6962SLisandro Dalcin start += count; 148492bb6962SLisandro Dalcin } 148592bb6962SLisandro Dalcin 148692bb6962SLisandro Dalcin } else { 148792bb6962SLisandro Dalcin 148892bb6962SLisandro Dalcin /* Partitioning by adjacency of diagonal block */ 148992bb6962SLisandro Dalcin 149092bb6962SLisandro Dalcin const PetscInt *numbering; 149192bb6962SLisandro Dalcin PetscInt *count,nidx,*indices,*newidx,start=0; 149292bb6962SLisandro Dalcin /* Get node count in each partition */ 14939566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n,&count)); 14949566063dSJacob Faibussowitsch PetscCall(ISPartitioningCount(ispart,n,count)); 149592bb6962SLisandro Dalcin if (isbaij && bs > 1) { /* adjust for the block-aij case */ 149692bb6962SLisandro Dalcin for (i=0; i<n; i++) count[i] *= bs; 149792bb6962SLisandro Dalcin } 149892bb6962SLisandro Dalcin /* Build indices from node numbering */ 14999566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(isnumb,&nidx)); 15009566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nidx,&indices)); 150192bb6962SLisandro Dalcin for (i=0; i<nidx; i++) indices[i] = i; /* needs to be initialized */ 15029566063dSJacob Faibussowitsch PetscCall(ISGetIndices(isnumb,&numbering)); 15039566063dSJacob Faibussowitsch PetscCall(PetscSortIntWithPermutation(nidx,numbering,indices)); 15049566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(isnumb,&numbering)); 150592bb6962SLisandro Dalcin if (isbaij && bs > 1) { /* adjust for the block-aij case */ 15069566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nidx*bs,&newidx)); 15072fa5cd67SKarl Rupp for (i=0; i<nidx; i++) { 15082fa5cd67SKarl Rupp for (j=0; j<bs; j++) newidx[i*bs+j] = indices[i]*bs + j; 15092fa5cd67SKarl Rupp } 15109566063dSJacob Faibussowitsch PetscCall(PetscFree(indices)); 151192bb6962SLisandro Dalcin nidx *= bs; 151292bb6962SLisandro Dalcin indices = newidx; 151392bb6962SLisandro Dalcin } 151492bb6962SLisandro Dalcin /* Shift to get global indices */ 151592bb6962SLisandro Dalcin for (i=0; i<nidx; i++) indices[i] += rstart; 151692bb6962SLisandro Dalcin 151792bb6962SLisandro Dalcin /* Build the index sets for each block */ 151892bb6962SLisandro Dalcin for (i=0; i<n; i++) { 15199566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF,count[i],&indices[start],PETSC_COPY_VALUES,&is[i])); 15209566063dSJacob Faibussowitsch PetscCall(ISSort(is[i])); 152192bb6962SLisandro Dalcin start += count[i]; 152292bb6962SLisandro Dalcin } 152392bb6962SLisandro Dalcin 15249566063dSJacob Faibussowitsch PetscCall(PetscFree(count)); 15259566063dSJacob Faibussowitsch PetscCall(PetscFree(indices)); 15269566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isnumb)); 15279566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ispart)); 152892bb6962SLisandro Dalcin 152992bb6962SLisandro Dalcin } 153092bb6962SLisandro Dalcin PetscFunctionReturn(0); 153192bb6962SLisandro Dalcin } 153292bb6962SLisandro Dalcin 153392bb6962SLisandro Dalcin /*@C 153492bb6962SLisandro Dalcin PCASMDestroySubdomains - Destroys the index sets created with 153592bb6962SLisandro Dalcin PCASMCreateSubdomains(). Should be called after setting subdomains 153692bb6962SLisandro Dalcin with PCASMSetLocalSubdomains(). 153792bb6962SLisandro Dalcin 153892bb6962SLisandro Dalcin Collective 153992bb6962SLisandro Dalcin 154092bb6962SLisandro Dalcin Input Parameters: 154192bb6962SLisandro Dalcin + n - the number of index sets 15422b691e39SMatthew Knepley . is - the array of index sets 15430298fd71SBarry Smith - is_local - the array of local index sets, can be NULL 154492bb6962SLisandro Dalcin 154592bb6962SLisandro Dalcin Level: advanced 154692bb6962SLisandro Dalcin 1547db781477SPatrick Sanan .seealso: `PCASMCreateSubdomains()`, `PCASMSetLocalSubdomains()` 154892bb6962SLisandro Dalcin @*/ 15497087cfbeSBarry Smith PetscErrorCode PCASMDestroySubdomains(PetscInt n, IS is[], IS is_local[]) 155092bb6962SLisandro Dalcin { 155192bb6962SLisandro Dalcin PetscInt i; 15525fd66863SKarl Rupp 155392bb6962SLisandro Dalcin PetscFunctionBegin; 1554a35b7b57SDmitry Karpeev if (n <= 0) PetscFunctionReturn(0); 1555a35b7b57SDmitry Karpeev if (is) { 155692bb6962SLisandro Dalcin PetscValidPointer(is,2); 15579566063dSJacob Faibussowitsch for (i=0; i<n; i++) PetscCall(ISDestroy(&is[i])); 15589566063dSJacob Faibussowitsch PetscCall(PetscFree(is)); 1559a35b7b57SDmitry Karpeev } 15602b691e39SMatthew Knepley if (is_local) { 15612b691e39SMatthew Knepley PetscValidPointer(is_local,3); 15629566063dSJacob Faibussowitsch for (i=0; i<n; i++) PetscCall(ISDestroy(&is_local[i])); 15639566063dSJacob Faibussowitsch PetscCall(PetscFree(is_local)); 15642b691e39SMatthew Knepley } 156592bb6962SLisandro Dalcin PetscFunctionReturn(0); 156692bb6962SLisandro Dalcin } 156792bb6962SLisandro Dalcin 15684b9ad928SBarry Smith /*@ 15694b9ad928SBarry Smith PCASMCreateSubdomains2D - Creates the index sets for the overlapping Schwarz 15704b9ad928SBarry Smith preconditioner for a two-dimensional problem on a regular grid. 15714b9ad928SBarry Smith 15724b9ad928SBarry Smith Not Collective 15734b9ad928SBarry Smith 15744b9ad928SBarry Smith Input Parameters: 15756b867d5aSJose E. Roman + m - the number of mesh points in the x direction 15766b867d5aSJose E. Roman . n - the number of mesh points in the y direction 15776b867d5aSJose E. Roman . M - the number of subdomains in the x direction 15786b867d5aSJose E. Roman . N - the number of subdomains in the y direction 15794b9ad928SBarry Smith . dof - degrees of freedom per node 15804b9ad928SBarry Smith - overlap - overlap in mesh lines 15814b9ad928SBarry Smith 15824b9ad928SBarry Smith Output Parameters: 15834b9ad928SBarry Smith + Nsub - the number of subdomains created 15843d03bb48SJed Brown . is - array of index sets defining overlapping (if overlap > 0) subdomains 15853d03bb48SJed Brown - is_local - array of index sets defining non-overlapping subdomains 15864b9ad928SBarry Smith 15874b9ad928SBarry Smith Note: 15884b9ad928SBarry Smith Presently PCAMSCreateSubdomains2d() is valid only for sequential 15894b9ad928SBarry Smith preconditioners. More general related routines are 15904b9ad928SBarry Smith PCASMSetTotalSubdomains() and PCASMSetLocalSubdomains(). 15914b9ad928SBarry Smith 15924b9ad928SBarry Smith Level: advanced 15934b9ad928SBarry Smith 1594db781477SPatrick Sanan .seealso: `PCASMSetTotalSubdomains()`, `PCASMSetLocalSubdomains()`, `PCASMGetSubKSP()`, 1595db781477SPatrick Sanan `PCASMSetOverlap()` 15964b9ad928SBarry Smith @*/ 15977087cfbeSBarry Smith PetscErrorCode PCASMCreateSubdomains2D(PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt dof,PetscInt overlap,PetscInt *Nsub,IS **is,IS **is_local) 15984b9ad928SBarry Smith { 15993d03bb48SJed Brown PetscInt i,j,height,width,ystart,xstart,yleft,yright,xleft,xright,loc_outer; 160013f74950SBarry Smith PetscInt nidx,*idx,loc,ii,jj,count; 16014b9ad928SBarry Smith 16024b9ad928SBarry Smith PetscFunctionBegin; 16037827d75bSBarry Smith PetscCheck(dof == 1,PETSC_COMM_SELF,PETSC_ERR_SUP,"dof must be 1"); 16044b9ad928SBarry Smith 16054b9ad928SBarry Smith *Nsub = N*M; 16069566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(*Nsub,is)); 16079566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(*Nsub,is_local)); 16084b9ad928SBarry Smith ystart = 0; 16093d03bb48SJed Brown loc_outer = 0; 16104b9ad928SBarry Smith for (i=0; i<N; i++) { 16114b9ad928SBarry Smith height = n/N + ((n % N) > i); /* height of subdomain */ 16127827d75bSBarry Smith PetscCheck(height >= 2,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many N subdomains for mesh dimension n"); 16134b9ad928SBarry Smith yleft = ystart - overlap; if (yleft < 0) yleft = 0; 16144b9ad928SBarry Smith yright = ystart + height + overlap; if (yright > n) yright = n; 16154b9ad928SBarry Smith xstart = 0; 16164b9ad928SBarry Smith for (j=0; j<M; j++) { 16174b9ad928SBarry Smith width = m/M + ((m % M) > j); /* width of subdomain */ 16187827d75bSBarry Smith PetscCheck(width >= 2,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many M subdomains for mesh dimension m"); 16194b9ad928SBarry Smith xleft = xstart - overlap; if (xleft < 0) xleft = 0; 16204b9ad928SBarry Smith xright = xstart + width + overlap; if (xright > m) xright = m; 16214b9ad928SBarry Smith nidx = (xright - xleft)*(yright - yleft); 16229566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nidx,&idx)); 16234b9ad928SBarry Smith loc = 0; 16244b9ad928SBarry Smith for (ii=yleft; ii<yright; ii++) { 16254b9ad928SBarry Smith count = m*ii + xleft; 16262fa5cd67SKarl Rupp for (jj=xleft; jj<xright; jj++) idx[loc++] = count++; 16274b9ad928SBarry Smith } 16289566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF,nidx,idx,PETSC_COPY_VALUES,(*is)+loc_outer)); 16293d03bb48SJed Brown if (overlap == 0) { 16309566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)(*is)[loc_outer])); 16312fa5cd67SKarl Rupp 16323d03bb48SJed Brown (*is_local)[loc_outer] = (*is)[loc_outer]; 16333d03bb48SJed Brown } else { 16343d03bb48SJed Brown for (loc=0,ii=ystart; ii<ystart+height; ii++) { 16353d03bb48SJed Brown for (jj=xstart; jj<xstart+width; jj++) { 16363d03bb48SJed Brown idx[loc++] = m*ii + jj; 16373d03bb48SJed Brown } 16383d03bb48SJed Brown } 16399566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF,loc,idx,PETSC_COPY_VALUES,*is_local+loc_outer)); 16403d03bb48SJed Brown } 16419566063dSJacob Faibussowitsch PetscCall(PetscFree(idx)); 16424b9ad928SBarry Smith xstart += width; 16433d03bb48SJed Brown loc_outer++; 16444b9ad928SBarry Smith } 16454b9ad928SBarry Smith ystart += height; 16464b9ad928SBarry Smith } 16479566063dSJacob Faibussowitsch for (i=0; i<*Nsub; i++) PetscCall(ISSort((*is)[i])); 16484b9ad928SBarry Smith PetscFunctionReturn(0); 16494b9ad928SBarry Smith } 16504b9ad928SBarry Smith 16514b9ad928SBarry Smith /*@C 16524b9ad928SBarry Smith PCASMGetLocalSubdomains - Gets the local subdomains (for this processor 16534b9ad928SBarry Smith only) for the additive Schwarz preconditioner. 16544b9ad928SBarry Smith 1655ad4df100SBarry Smith Not Collective 16564b9ad928SBarry Smith 16574b9ad928SBarry Smith Input Parameter: 16584b9ad928SBarry Smith . pc - the preconditioner context 16594b9ad928SBarry Smith 16604b9ad928SBarry Smith Output Parameters: 1661f17a6784SPierre Jolivet + n - if requested, the number of subdomains for this processor (default value = 1) 1662f17a6784SPierre Jolivet . is - if requested, the index sets that define the subdomains for this processor 1663f17a6784SPierre Jolivet - is_local - if requested, the index sets that define the local part of the subdomains for this processor (can be NULL) 16644b9ad928SBarry Smith 16654b9ad928SBarry Smith Notes: 16664b9ad928SBarry Smith The IS numbering is in the parallel, global numbering of the vector. 16674b9ad928SBarry Smith 16684b9ad928SBarry Smith Level: advanced 16694b9ad928SBarry Smith 1670db781477SPatrick Sanan .seealso: `PCASMSetTotalSubdomains()`, `PCASMSetOverlap()`, `PCASMGetSubKSP()`, 1671db781477SPatrick Sanan `PCASMCreateSubdomains2D()`, `PCASMSetLocalSubdomains()`, `PCASMGetLocalSubmatrices()` 16724b9ad928SBarry Smith @*/ 16737087cfbeSBarry Smith PetscErrorCode PCASMGetLocalSubdomains(PC pc,PetscInt *n,IS *is[],IS *is_local[]) 16744b9ad928SBarry Smith { 16752a808120SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 1676ace3abfcSBarry Smith PetscBool match; 16774b9ad928SBarry Smith 16784b9ad928SBarry Smith PetscFunctionBegin; 16790700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1680f17a6784SPierre Jolivet if (n) PetscValidIntPointer(n,2); 168192bb6962SLisandro Dalcin if (is) PetscValidPointer(is,3); 1682f17a6784SPierre Jolivet if (is_local) PetscValidPointer(is_local,4); 16839566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc,PCASM,&match)); 168428b400f6SJacob Faibussowitsch PetscCheck(match,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"PC is not a PCASM"); 16854b9ad928SBarry Smith if (n) *n = osm->n_local_true; 16864b9ad928SBarry Smith if (is) *is = osm->is; 16872b691e39SMatthew Knepley if (is_local) *is_local = osm->is_local; 16884b9ad928SBarry Smith PetscFunctionReturn(0); 16894b9ad928SBarry Smith } 16904b9ad928SBarry Smith 16914b9ad928SBarry Smith /*@C 16924b9ad928SBarry Smith PCASMGetLocalSubmatrices - Gets the local submatrices (for this processor 16934b9ad928SBarry Smith only) for the additive Schwarz preconditioner. 16944b9ad928SBarry Smith 1695ad4df100SBarry Smith Not Collective 16964b9ad928SBarry Smith 16974b9ad928SBarry Smith Input Parameter: 16984b9ad928SBarry Smith . pc - the preconditioner context 16994b9ad928SBarry Smith 17004b9ad928SBarry Smith Output Parameters: 1701f17a6784SPierre Jolivet + n - if requested, the number of matrices for this processor (default value = 1) 1702f17a6784SPierre Jolivet - mat - if requested, the matrices 17034b9ad928SBarry Smith 17044b9ad928SBarry Smith Level: advanced 17054b9ad928SBarry Smith 170695452b02SPatrick Sanan Notes: 170734a84908Sprj- Call after PCSetUp() (or KSPSetUp()) but before PCApply() and before PCSetUpOnBlocks()) 1708cf739d55SBarry Smith 170934a84908Sprj- Usually one would use PCSetModifySubMatrices() to change the submatrices in building the preconditioner. 1710cf739d55SBarry Smith 1711db781477SPatrick Sanan .seealso: `PCASMSetTotalSubdomains()`, `PCASMSetOverlap()`, `PCASMGetSubKSP()`, 1712db781477SPatrick Sanan `PCASMCreateSubdomains2D()`, `PCASMSetLocalSubdomains()`, `PCASMGetLocalSubdomains()`, `PCSetModifySubMatrices()` 17134b9ad928SBarry Smith @*/ 17147087cfbeSBarry Smith PetscErrorCode PCASMGetLocalSubmatrices(PC pc,PetscInt *n,Mat *mat[]) 17154b9ad928SBarry Smith { 17164b9ad928SBarry Smith PC_ASM *osm; 1717ace3abfcSBarry Smith PetscBool match; 17184b9ad928SBarry Smith 17194b9ad928SBarry Smith PetscFunctionBegin; 17200700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1721f17a6784SPierre Jolivet if (n) PetscValidIntPointer(n,2); 172292bb6962SLisandro Dalcin if (mat) PetscValidPointer(mat,3); 172328b400f6SJacob Faibussowitsch PetscCheck(pc->setupcalled,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call after KSPSetUp() or PCSetUp()."); 17249566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc,PCASM,&match)); 172592bb6962SLisandro Dalcin if (!match) { 172692bb6962SLisandro Dalcin if (n) *n = 0; 17270298fd71SBarry Smith if (mat) *mat = NULL; 172892bb6962SLisandro Dalcin } else { 17294b9ad928SBarry Smith osm = (PC_ASM*)pc->data; 17304b9ad928SBarry Smith if (n) *n = osm->n_local_true; 17314b9ad928SBarry Smith if (mat) *mat = osm->pmat; 173292bb6962SLisandro Dalcin } 17334b9ad928SBarry Smith PetscFunctionReturn(0); 17344b9ad928SBarry Smith } 1735d709ea83SDmitry Karpeev 1736d709ea83SDmitry Karpeev /*@ 1737d709ea83SDmitry Karpeev PCASMSetDMSubdomains - Indicates whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible. 1738f1ee410cSBarry Smith 1739d709ea83SDmitry Karpeev Logically Collective 1740d709ea83SDmitry Karpeev 1741d8d19677SJose E. Roman Input Parameters: 1742d709ea83SDmitry Karpeev + pc - the preconditioner 1743d709ea83SDmitry Karpeev - flg - boolean indicating whether to use subdomains defined by the DM 1744d709ea83SDmitry Karpeev 1745d709ea83SDmitry Karpeev Options Database Key: 1746147403d9SBarry Smith . -pc_asm_dm_subdomains <bool> - use subdomains defined by the DM 1747d709ea83SDmitry Karpeev 1748d709ea83SDmitry Karpeev Level: intermediate 1749d709ea83SDmitry Karpeev 1750d709ea83SDmitry Karpeev Notes: 1751d709ea83SDmitry Karpeev PCASMSetTotalSubdomains() and PCASMSetOverlap() take precedence over PCASMSetDMSubdomains(), 1752d709ea83SDmitry Karpeev so setting either of the first two effectively turns the latter off. 1753d709ea83SDmitry Karpeev 1754db781477SPatrick Sanan .seealso: `PCASMGetDMSubdomains()`, `PCASMSetTotalSubdomains()`, `PCASMSetOverlap()` 1755db781477SPatrick Sanan `PCASMCreateSubdomains2D()`, `PCASMSetLocalSubdomains()`, `PCASMGetLocalSubdomains()` 1756d709ea83SDmitry Karpeev @*/ 1757d709ea83SDmitry Karpeev PetscErrorCode PCASMSetDMSubdomains(PC pc,PetscBool flg) 1758d709ea83SDmitry Karpeev { 1759d709ea83SDmitry Karpeev PC_ASM *osm = (PC_ASM*)pc->data; 1760d709ea83SDmitry Karpeev PetscBool match; 1761d709ea83SDmitry Karpeev 1762d709ea83SDmitry Karpeev PetscFunctionBegin; 1763d709ea83SDmitry Karpeev PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1764d709ea83SDmitry Karpeev PetscValidLogicalCollectiveBool(pc,flg,2); 176528b400f6SJacob Faibussowitsch PetscCheck(!pc->setupcalled,((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Not for a setup PC."); 17669566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc,PCASM,&match)); 1767d709ea83SDmitry Karpeev if (match) { 1768d709ea83SDmitry Karpeev osm->dm_subdomains = flg; 1769d709ea83SDmitry Karpeev } 1770d709ea83SDmitry Karpeev PetscFunctionReturn(0); 1771d709ea83SDmitry Karpeev } 1772d709ea83SDmitry Karpeev 1773d709ea83SDmitry Karpeev /*@ 1774d709ea83SDmitry Karpeev PCASMGetDMSubdomains - Returns flag indicating whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible. 1775d709ea83SDmitry Karpeev Not Collective 1776d709ea83SDmitry Karpeev 1777d709ea83SDmitry Karpeev Input Parameter: 1778d709ea83SDmitry Karpeev . pc - the preconditioner 1779d709ea83SDmitry Karpeev 1780d709ea83SDmitry Karpeev Output Parameter: 1781d709ea83SDmitry Karpeev . flg - boolean indicating whether to use subdomains defined by the DM 1782d709ea83SDmitry Karpeev 1783d709ea83SDmitry Karpeev Level: intermediate 1784d709ea83SDmitry Karpeev 1785db781477SPatrick Sanan .seealso: `PCASMSetDMSubdomains()`, `PCASMSetTotalSubdomains()`, `PCASMSetOverlap()` 1786db781477SPatrick Sanan `PCASMCreateSubdomains2D()`, `PCASMSetLocalSubdomains()`, `PCASMGetLocalSubdomains()` 1787d709ea83SDmitry Karpeev @*/ 1788d709ea83SDmitry Karpeev PetscErrorCode PCASMGetDMSubdomains(PC pc,PetscBool* flg) 1789d709ea83SDmitry Karpeev { 1790d709ea83SDmitry Karpeev PC_ASM *osm = (PC_ASM*)pc->data; 1791d709ea83SDmitry Karpeev PetscBool match; 1792d709ea83SDmitry Karpeev 1793d709ea83SDmitry Karpeev PetscFunctionBegin; 1794d709ea83SDmitry Karpeev PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1795534a8f05SLisandro Dalcin PetscValidBoolPointer(flg,2); 17969566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc,PCASM,&match)); 179756ea09ceSStefano Zampini if (match) *flg = osm->dm_subdomains; 179856ea09ceSStefano Zampini else *flg = PETSC_FALSE; 1799d709ea83SDmitry Karpeev PetscFunctionReturn(0); 1800d709ea83SDmitry Karpeev } 180180ec0b7dSPatrick Sanan 180280ec0b7dSPatrick Sanan /*@ 180380ec0b7dSPatrick Sanan PCASMGetSubMatType - Gets the matrix type used for ASM subsolves, as a string. 180480ec0b7dSPatrick Sanan 180580ec0b7dSPatrick Sanan Not Collective 180680ec0b7dSPatrick Sanan 180780ec0b7dSPatrick Sanan Input Parameter: 180880ec0b7dSPatrick Sanan . pc - the PC 180980ec0b7dSPatrick Sanan 181080ec0b7dSPatrick Sanan Output Parameter: 1811f1ee410cSBarry Smith . -pc_asm_sub_mat_type - name of matrix type 181280ec0b7dSPatrick Sanan 181380ec0b7dSPatrick Sanan Level: advanced 181480ec0b7dSPatrick Sanan 1815db781477SPatrick Sanan .seealso: `PCASMSetSubMatType()`, `PCASM`, `PCSetType()`, `VecSetType()`, `MatType`, `Mat` 181680ec0b7dSPatrick Sanan @*/ 181756ea09ceSStefano Zampini PetscErrorCode PCASMGetSubMatType(PC pc,MatType *sub_mat_type) 181856ea09ceSStefano Zampini { 181956ea09ceSStefano Zampini PetscFunctionBegin; 182056ea09ceSStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1821cac4c232SBarry Smith PetscTryMethod(pc,"PCASMGetSubMatType_C",(PC,MatType*),(pc,sub_mat_type)); 182280ec0b7dSPatrick Sanan PetscFunctionReturn(0); 182380ec0b7dSPatrick Sanan } 182480ec0b7dSPatrick Sanan 182580ec0b7dSPatrick Sanan /*@ 182680ec0b7dSPatrick Sanan PCASMSetSubMatType - Set the type of matrix used for ASM subsolves 182780ec0b7dSPatrick Sanan 182880ec0b7dSPatrick Sanan Collective on Mat 182980ec0b7dSPatrick Sanan 183080ec0b7dSPatrick Sanan Input Parameters: 183180ec0b7dSPatrick Sanan + pc - the PC object 183280ec0b7dSPatrick Sanan - sub_mat_type - matrix type 183380ec0b7dSPatrick Sanan 183480ec0b7dSPatrick Sanan Options Database Key: 183580ec0b7dSPatrick Sanan . -pc_asm_sub_mat_type <sub_mat_type> - Sets the matrix type used for subsolves, for example, seqaijviennacl. If you specify a base name like aijviennacl, the corresponding sequential type is assumed. 183680ec0b7dSPatrick Sanan 183780ec0b7dSPatrick Sanan Notes: 183880ec0b7dSPatrick Sanan See "${PETSC_DIR}/include/petscmat.h" for available types 183980ec0b7dSPatrick Sanan 184080ec0b7dSPatrick Sanan Level: advanced 184180ec0b7dSPatrick Sanan 1842db781477SPatrick Sanan .seealso: `PCASMGetSubMatType()`, `PCASM`, `PCSetType()`, `VecSetType()`, `MatType`, `Mat` 184380ec0b7dSPatrick Sanan @*/ 184480ec0b7dSPatrick Sanan PetscErrorCode PCASMSetSubMatType(PC pc,MatType sub_mat_type) 184580ec0b7dSPatrick Sanan { 184656ea09ceSStefano Zampini PetscFunctionBegin; 184756ea09ceSStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1848cac4c232SBarry Smith PetscTryMethod(pc,"PCASMSetSubMatType_C",(PC,MatType),(pc,sub_mat_type)); 184980ec0b7dSPatrick Sanan PetscFunctionReturn(0); 185080ec0b7dSPatrick Sanan } 1851