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*/ 144b9ad928SBarry Smith 156849ba73SBarry Smith static PetscErrorCode PCView_ASM(PC pc,PetscViewer viewer) 164b9ad928SBarry Smith { 1792bb6962SLisandro Dalcin PC_ASM *osm = (PC_ASM*)pc->data; 1813f74950SBarry Smith PetscMPIInt rank; 194d219a6aSLisandro Dalcin PetscInt i,bsz; 20ace3abfcSBarry Smith PetscBool iascii,isstring; 214b9ad928SBarry Smith PetscViewer sviewer; 22ed155784SPierre Jolivet PetscViewerFormat format; 23ed155784SPierre Jolivet const char *prefix; 244b9ad928SBarry Smith 254b9ad928SBarry Smith PetscFunctionBegin; 269566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii)); 279566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring)); 2832077d6dSBarry Smith if (iascii) { 293d03bb48SJed Brown char overlaps[256] = "user-defined overlap",blocks[256] = "total subdomain blocks not yet set"; 3063a3b9bcSJacob Faibussowitsch if (osm->overlap >= 0) PetscCall(PetscSNPrintf(overlaps,sizeof(overlaps),"amount of overlap = %" PetscInt_FMT,osm->overlap)); 3163a3b9bcSJacob Faibussowitsch if (osm->n > 0) PetscCall(PetscSNPrintf(blocks,sizeof(blocks),"total subdomain blocks = %" PetscInt_FMT,osm->n)); 329566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," %s, %s\n",blocks,overlaps)); 339566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," restriction/interpolation type - %s\n",PCASMTypes[osm->type])); 349566063dSJacob Faibussowitsch if (osm->dm_subdomains) PetscCall(PetscViewerASCIIPrintf(viewer," Additive Schwarz: using DM to define subdomains\n")); 359566063dSJacob Faibussowitsch if (osm->loctype != PC_COMPOSITE_ADDITIVE) PetscCall(PetscViewerASCIIPrintf(viewer," Additive Schwarz: local solve composition type - %s\n",PCCompositeTypes[osm->loctype])); 369566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank)); 379566063dSJacob Faibussowitsch PetscCall(PetscViewerGetFormat(viewer,&format)); 38ed155784SPierre Jolivet if (format != PETSC_VIEWER_ASCII_INFO_DETAIL) { 394d219a6aSLisandro Dalcin if (osm->ksp) { 409566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Local solver information for first block is in the following KSP and PC objects on rank 0:\n")); 419566063dSJacob Faibussowitsch PetscCall(PCGetOptionsPrefix(pc,&prefix)); 429566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Use -%sksp_view ::ascii_info_detail to display information for all blocks\n",prefix?prefix:"")); 439566063dSJacob Faibussowitsch PetscCall(PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer)); 44dd400576SPatrick Sanan if (rank == 0) { 459566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 469566063dSJacob Faibussowitsch PetscCall(KSPView(osm->ksp[0],sviewer)); 479566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 484b9ad928SBarry Smith } 499566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer)); 504d219a6aSLisandro Dalcin } 514b9ad928SBarry Smith } else { 529566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 5363a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(viewer," [%d] number of local blocks = %" PetscInt_FMT "\n",(int)rank,osm->n_local_true)); 549566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 559566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Local solver information for each block is in the following KSP and PC objects:\n")); 569566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 579566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"- - - - - - - - - - - - - - - - - -\n")); 589566063dSJacob Faibussowitsch PetscCall(PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer)); 59dd2fa690SBarry Smith for (i=0; i<osm->n_local_true; i++) { 609566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(osm->is[i],&bsz)); 6163a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(sviewer,"[%d] local block number %" PetscInt_FMT ", size = %" PetscInt_FMT "\n",(int)rank,i,bsz)); 629566063dSJacob Faibussowitsch PetscCall(KSPView(osm->ksp[i],sviewer)); 639566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(sviewer,"- - - - - - - - - - - - - - - - - -\n")); 644b9ad928SBarry Smith } 659566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer)); 669566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 679566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 689566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 694b9ad928SBarry Smith } 704b9ad928SBarry Smith } else if (isstring) { 7163a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerStringSPrintf(viewer," blocks=%" PetscInt_FMT ", overlap=%" PetscInt_FMT ", type=%s",osm->n,osm->overlap,PCASMTypes[osm->type])); 729566063dSJacob Faibussowitsch PetscCall(PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer)); 739566063dSJacob Faibussowitsch if (osm->ksp) PetscCall(KSPView(osm->ksp[0],sviewer)); 749566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer)); 754b9ad928SBarry Smith } 764b9ad928SBarry Smith PetscFunctionReturn(0); 774b9ad928SBarry Smith } 784b9ad928SBarry Smith 7992bb6962SLisandro Dalcin static PetscErrorCode PCASMPrintSubdomains(PC pc) 8092bb6962SLisandro Dalcin { 8192bb6962SLisandro Dalcin PC_ASM *osm = (PC_ASM*)pc->data; 8292bb6962SLisandro Dalcin const char *prefix; 8392bb6962SLisandro Dalcin char fname[PETSC_MAX_PATH_LEN+1]; 84643535aeSDmitry Karpeev PetscViewer viewer, sviewer; 85643535aeSDmitry Karpeev char *s; 8692bb6962SLisandro Dalcin PetscInt i,j,nidx; 8792bb6962SLisandro Dalcin const PetscInt *idx; 88643535aeSDmitry Karpeev PetscMPIInt rank, size; 89846783a0SBarry Smith 9092bb6962SLisandro Dalcin PetscFunctionBegin; 919566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size)); 929566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank)); 939566063dSJacob Faibussowitsch PetscCall(PCGetOptionsPrefix(pc,&prefix)); 949566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL,prefix,"-pc_asm_print_subdomains",fname,sizeof(fname),NULL)); 959566063dSJacob Faibussowitsch if (fname[0] == 0) PetscCall(PetscStrcpy(fname,"stdout")); 969566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIOpen(PetscObjectComm((PetscObject)pc),fname,&viewer)); 97643535aeSDmitry Karpeev for (i=0; i<osm->n_local; i++) { 98643535aeSDmitry Karpeev if (i < osm->n_local_true) { 999566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(osm->is[i],&nidx)); 1009566063dSJacob Faibussowitsch PetscCall(ISGetIndices(osm->is[i],&idx)); 101643535aeSDmitry Karpeev /* Print to a string viewer; no more than 15 characters per index plus 512 char for the header.*/ 10236a9e3b9SBarry Smith #define len 16*(nidx+1)+512 1039566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(len,&s)); 1049566063dSJacob Faibussowitsch PetscCall(PetscViewerStringOpen(PETSC_COMM_SELF, s, len, &sviewer)); 10536a9e3b9SBarry Smith #undef len 10663a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerStringSPrintf(sviewer, "[%d:%d] Subdomain %" PetscInt_FMT " with overlap:\n", rank, size, i)); 10792bb6962SLisandro Dalcin for (j=0; j<nidx; j++) { 10863a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerStringSPrintf(sviewer,"%" PetscInt_FMT " ",idx[j])); 10992bb6962SLisandro Dalcin } 1109566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(osm->is[i],&idx)); 1119566063dSJacob Faibussowitsch PetscCall(PetscViewerStringSPrintf(sviewer,"\n")); 1129566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&sviewer)); 1139566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 11463a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "%s", s)); 1159566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 1169566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 1179566063dSJacob Faibussowitsch PetscCall(PetscFree(s)); 1182b691e39SMatthew Knepley if (osm->is_local) { 119643535aeSDmitry Karpeev /* Print to a string viewer; no more than 15 characters per index plus 512 char for the header.*/ 12036a9e3b9SBarry Smith #define len 16*(nidx+1)+512 1219566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(len, &s)); 1229566063dSJacob Faibussowitsch PetscCall(PetscViewerStringOpen(PETSC_COMM_SELF, s, len, &sviewer)); 12336a9e3b9SBarry Smith #undef len 12463a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerStringSPrintf(sviewer, "[%d:%d] Subdomain %" PetscInt_FMT " without overlap:\n", rank, size, i)); 1259566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(osm->is_local[i],&nidx)); 1269566063dSJacob Faibussowitsch PetscCall(ISGetIndices(osm->is_local[i],&idx)); 1272b691e39SMatthew Knepley for (j=0; j<nidx; j++) { 12863a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerStringSPrintf(sviewer,"%" PetscInt_FMT " ",idx[j])); 1292b691e39SMatthew Knepley } 1309566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(osm->is_local[i],&idx)); 1319566063dSJacob Faibussowitsch PetscCall(PetscViewerStringSPrintf(sviewer,"\n")); 1329566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&sviewer)); 1339566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 13463a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "%s", s)); 1359566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 1369566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 1379566063dSJacob Faibussowitsch PetscCall(PetscFree(s)); 138643535aeSDmitry Karpeev } 1392fa5cd67SKarl Rupp } else { 140643535aeSDmitry Karpeev /* Participate in collective viewer calls. */ 1419566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 1429566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 1439566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 144643535aeSDmitry Karpeev /* Assume either all ranks have is_local or none do. */ 145643535aeSDmitry Karpeev if (osm->is_local) { 1469566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 1479566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 1489566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 149643535aeSDmitry Karpeev } 150fdc48646SDmitry Karpeev } 15192bb6962SLisandro Dalcin } 1529566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 1539566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 15492bb6962SLisandro Dalcin PetscFunctionReturn(0); 15592bb6962SLisandro Dalcin } 15692bb6962SLisandro Dalcin 1576849ba73SBarry Smith static PetscErrorCode PCSetUp_ASM(PC pc) 1584b9ad928SBarry Smith { 1594b9ad928SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 16057501b6eSBarry Smith PetscBool flg; 16187e86712SBarry Smith PetscInt i,m,m_local; 1624b9ad928SBarry Smith MatReuse scall = MAT_REUSE_MATRIX; 1634b9ad928SBarry Smith IS isl; 1644b9ad928SBarry Smith KSP ksp; 1654b9ad928SBarry Smith PC subpc; 1662dcb1b2aSMatthew Knepley const char *prefix,*pprefix; 16723ce1328SBarry Smith Vec vec; 1680298fd71SBarry Smith DM *domain_dm = NULL; 1694b9ad928SBarry Smith 1704b9ad928SBarry Smith PetscFunctionBegin; 1714b9ad928SBarry Smith if (!pc->setupcalled) { 172265a2120SBarry Smith PetscInt m; 17392bb6962SLisandro Dalcin 1742b837212SDmitry Karpeev /* Note: if subdomains have been set either via PCASMSetTotalSubdomains() or via PCASMSetLocalSubdomains(), osm->n_local_true will not be PETSC_DECIDE */ 1752b837212SDmitry Karpeev if (osm->n_local_true == PETSC_DECIDE) { 17669ca1f37SDmitry Karpeev /* no subdomains given */ 17765db9045SDmitry Karpeev /* try pc->dm first, if allowed */ 178d709ea83SDmitry Karpeev if (osm->dm_subdomains && pc->dm) { 179feb221f8SDmitry Karpeev PetscInt num_domains, d; 180feb221f8SDmitry Karpeev char **domain_names; 1818d4ac253SDmitry Karpeev IS *inner_domain_is, *outer_domain_is; 1829566063dSJacob Faibussowitsch PetscCall(DMCreateDomainDecomposition(pc->dm, &num_domains, &domain_names, &inner_domain_is, &outer_domain_is, &domain_dm)); 183704f0395SPatrick Sanan osm->overlap = -1; /* We do not want to increase the overlap of the IS. 184704f0395SPatrick Sanan A future improvement of this code might allow one to use 185704f0395SPatrick Sanan DM-defined subdomains and also increase the overlap, 186704f0395SPatrick Sanan but that is not currently supported */ 18769ca1f37SDmitry Karpeev if (num_domains) { 1889566063dSJacob Faibussowitsch PetscCall(PCASMSetLocalSubdomains(pc, num_domains, outer_domain_is, inner_domain_is)); 18969ca1f37SDmitry Karpeev } 190feb221f8SDmitry Karpeev for (d = 0; d < num_domains; ++d) { 1919566063dSJacob Faibussowitsch if (domain_names) PetscCall(PetscFree(domain_names[d])); 1929566063dSJacob Faibussowitsch if (inner_domain_is) PetscCall(ISDestroy(&inner_domain_is[d])); 1939566063dSJacob Faibussowitsch if (outer_domain_is) PetscCall(ISDestroy(&outer_domain_is[d])); 194feb221f8SDmitry Karpeev } 1959566063dSJacob Faibussowitsch PetscCall(PetscFree(domain_names)); 1969566063dSJacob Faibussowitsch PetscCall(PetscFree(inner_domain_is)); 1979566063dSJacob Faibussowitsch PetscCall(PetscFree(outer_domain_is)); 198feb221f8SDmitry Karpeev } 1992b837212SDmitry Karpeev if (osm->n_local_true == PETSC_DECIDE) { 20069ca1f37SDmitry Karpeev /* still no subdomains; use one subdomain per processor */ 2012b837212SDmitry Karpeev osm->n_local_true = 1; 202feb221f8SDmitry Karpeev } 2032b837212SDmitry Karpeev } 2042b837212SDmitry Karpeev { /* determine the global and max number of subdomains */ 2056ac3741eSJed Brown struct {PetscInt max,sum;} inwork,outwork; 206c10200c1SHong Zhang PetscMPIInt size; 207c10200c1SHong Zhang 2086ac3741eSJed Brown inwork.max = osm->n_local_true; 2096ac3741eSJed Brown inwork.sum = osm->n_local_true; 2101c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(&inwork,&outwork,1,MPIU_2INT,MPIU_MAXSUM_OP,PetscObjectComm((PetscObject)pc))); 2116ac3741eSJed Brown osm->n_local = outwork.max; 2126ac3741eSJed Brown osm->n = outwork.sum; 213c10200c1SHong Zhang 2149566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size)); 215c10200c1SHong Zhang if (outwork.max == 1 && outwork.sum == size) { 2167dae84e0SHong Zhang /* osm->n_local_true = 1 on all processes, set this option may enable use of optimized MatCreateSubMatrices() implementation */ 2179566063dSJacob Faibussowitsch PetscCall(MatSetOption(pc->pmat,MAT_SUBMAT_SINGLEIS,PETSC_TRUE)); 218c10200c1SHong Zhang } 2194b9ad928SBarry Smith } 22078904715SLisandro Dalcin if (!osm->is) { /* create the index sets */ 2219566063dSJacob Faibussowitsch PetscCall(PCASMCreateSubdomains(pc->pmat,osm->n_local_true,&osm->is)); 2224b9ad928SBarry Smith } 223f5234e35SJed Brown if (osm->n_local_true > 1 && !osm->is_local) { 2249566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(osm->n_local_true,&osm->is_local)); 225f5234e35SJed Brown for (i=0; i<osm->n_local_true; i++) { 226f5234e35SJed Brown if (osm->overlap > 0) { /* With positive overlap, osm->is[i] will be modified */ 2279566063dSJacob Faibussowitsch PetscCall(ISDuplicate(osm->is[i],&osm->is_local[i])); 2289566063dSJacob Faibussowitsch PetscCall(ISCopy(osm->is[i],osm->is_local[i])); 229f5234e35SJed Brown } else { 2309566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)osm->is[i])); 231f5234e35SJed Brown osm->is_local[i] = osm->is[i]; 232f5234e35SJed Brown } 233f5234e35SJed Brown } 234f5234e35SJed Brown } 2359566063dSJacob Faibussowitsch PetscCall(PCGetOptionsPrefix(pc,&prefix)); 2363d03bb48SJed Brown if (osm->overlap > 0) { 2374b9ad928SBarry Smith /* Extend the "overlapping" regions by a number of steps */ 2389566063dSJacob Faibussowitsch PetscCall(MatIncreaseOverlap(pc->pmat,osm->n_local_true,osm->is,osm->overlap)); 2393d03bb48SJed Brown } 2406ed231c7SMatthew Knepley if (osm->sort_indices) { 24192bb6962SLisandro Dalcin for (i=0; i<osm->n_local_true; i++) { 2429566063dSJacob Faibussowitsch PetscCall(ISSort(osm->is[i])); 2432b691e39SMatthew Knepley if (osm->is_local) { 2449566063dSJacob Faibussowitsch PetscCall(ISSort(osm->is_local[i])); 2452b691e39SMatthew Knepley } 2464b9ad928SBarry Smith } 2476ed231c7SMatthew Knepley } 24898d3e228SPierre Jolivet flg = PETSC_FALSE; 24998d3e228SPierre Jolivet PetscCall(PetscOptionsHasName(NULL,prefix,"-pc_asm_print_subdomains",&flg)); 25098d3e228SPierre Jolivet if (flg) PetscCall(PCASMPrintSubdomains(pc)); 251f6991133SBarry Smith if (!osm->ksp) { 25278904715SLisandro Dalcin /* Create the local solvers */ 2539566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(osm->n_local_true,&osm->ksp)); 254feb221f8SDmitry Karpeev if (domain_dm) { 2559566063dSJacob Faibussowitsch PetscCall(PetscInfo(pc,"Setting up ASM subproblems using the embedded DM\n")); 256feb221f8SDmitry Karpeev } 25792bb6962SLisandro Dalcin for (i=0; i<osm->n_local_true; i++) { 2589566063dSJacob Faibussowitsch PetscCall(KSPCreate(PETSC_COMM_SELF,&ksp)); 2599566063dSJacob Faibussowitsch PetscCall(KSPSetErrorIfNotConverged(ksp,pc->erroriffailure)); 2609566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)ksp)); 2619566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,1)); 2629566063dSJacob Faibussowitsch PetscCall(KSPSetType(ksp,KSPPREONLY)); 2639566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp,&subpc)); 2649566063dSJacob Faibussowitsch PetscCall(PCGetOptionsPrefix(pc,&prefix)); 2659566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(ksp,prefix)); 2669566063dSJacob Faibussowitsch PetscCall(KSPAppendOptionsPrefix(ksp,"sub_")); 267feb221f8SDmitry Karpeev if (domain_dm) { 2689566063dSJacob Faibussowitsch PetscCall(KSPSetDM(ksp, domain_dm[i])); 2699566063dSJacob Faibussowitsch PetscCall(KSPSetDMActive(ksp, PETSC_FALSE)); 2709566063dSJacob Faibussowitsch PetscCall(DMDestroy(&domain_dm[i])); 271feb221f8SDmitry Karpeev } 2724b9ad928SBarry Smith osm->ksp[i] = ksp; 2734b9ad928SBarry Smith } 274feb221f8SDmitry Karpeev if (domain_dm) { 2759566063dSJacob Faibussowitsch PetscCall(PetscFree(domain_dm)); 276feb221f8SDmitry Karpeev } 277f6991133SBarry Smith } 2781dd8081eSeaulisa 2799566063dSJacob Faibussowitsch PetscCall(ISConcatenate(PETSC_COMM_SELF, osm->n_local_true, osm->is, &osm->lis)); 2809566063dSJacob Faibussowitsch PetscCall(ISSortRemoveDups(osm->lis)); 2819566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(osm->lis, &m)); 2821dd8081eSeaulisa 2834b9ad928SBarry Smith scall = MAT_INITIAL_MATRIX; 2844b9ad928SBarry Smith } else { 2854b9ad928SBarry Smith /* 2864b9ad928SBarry Smith Destroy the blocks from the previous iteration 2874b9ad928SBarry Smith */ 288e09e08ccSBarry Smith if (pc->flag == DIFFERENT_NONZERO_PATTERN) { 2899566063dSJacob Faibussowitsch PetscCall(MatDestroyMatrices(osm->n_local_true,&osm->pmat)); 2904b9ad928SBarry Smith scall = MAT_INITIAL_MATRIX; 2914b9ad928SBarry Smith } 2924b9ad928SBarry Smith } 2934b9ad928SBarry Smith 294b58cb649SBarry Smith /* Destroy previous submatrices of a different type than pc->pmat since MAT_REUSE_MATRIX won't work in that case */ 295*82b5ce2aSStefano Zampini if (scall == MAT_REUSE_MATRIX && osm->sub_mat_type) { 296b58cb649SBarry Smith if (osm->n_local_true > 0) { 2979566063dSJacob Faibussowitsch PetscCall(MatDestroySubMatrices(osm->n_local_true,&osm->pmat)); 298b58cb649SBarry Smith } 299b58cb649SBarry Smith scall = MAT_INITIAL_MATRIX; 300b58cb649SBarry Smith } 301b58cb649SBarry Smith 30292bb6962SLisandro Dalcin /* 30392bb6962SLisandro Dalcin Extract out the submatrices 30492bb6962SLisandro Dalcin */ 3059566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(pc->pmat,osm->n_local_true,osm->is,osm->is,scall,&osm->pmat)); 30692bb6962SLisandro Dalcin if (scall == MAT_INITIAL_MATRIX) { 3079566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptionsPrefix((PetscObject)pc->pmat,&pprefix)); 30892bb6962SLisandro Dalcin for (i=0; i<osm->n_local_true; i++) { 3099566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)osm->pmat[i])); 3109566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)osm->pmat[i],pprefix)); 31192bb6962SLisandro Dalcin } 31292bb6962SLisandro Dalcin } 31380ec0b7dSPatrick Sanan 31480ec0b7dSPatrick Sanan /* Convert the types of the submatrices (if needbe) */ 31580ec0b7dSPatrick Sanan if (osm->sub_mat_type) { 31680ec0b7dSPatrick Sanan for (i=0; i<osm->n_local_true; i++) { 3179566063dSJacob Faibussowitsch PetscCall(MatConvert(osm->pmat[i],osm->sub_mat_type,MAT_INPLACE_MATRIX,&(osm->pmat[i]))); 31880ec0b7dSPatrick Sanan } 31980ec0b7dSPatrick Sanan } 32080ec0b7dSPatrick Sanan 32180ec0b7dSPatrick Sanan if (!pc->setupcalled) { 32256ea09ceSStefano Zampini VecType vtype; 32356ea09ceSStefano Zampini 32480ec0b7dSPatrick Sanan /* Create the local work vectors (from the local matrices) and scatter contexts */ 3259566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(pc->pmat,&vec,NULL)); 3265b3c0d42Seaulisa 3277827d75bSBarry 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()"); 3281dd8081eSeaulisa if (osm->is_local && osm->type == PC_ASM_RESTRICT && osm->loctype == PC_COMPOSITE_ADDITIVE) { 3299566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(osm->n_local_true,&osm->lprolongation)); 3301dd8081eSeaulisa } 3319566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(osm->n_local_true,&osm->lrestriction)); 3329566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(osm->n_local_true,&osm->x)); 3339566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(osm->n_local_true,&osm->y)); 334347574c9Seaulisa 3359566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(osm->lis,&m)); 3369566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF,m,0,1,&isl)); 3379566063dSJacob Faibussowitsch PetscCall(MatGetVecType(osm->pmat[0],&vtype)); 3389566063dSJacob Faibussowitsch PetscCall(VecCreate(PETSC_COMM_SELF,&osm->lx)); 3399566063dSJacob Faibussowitsch PetscCall(VecSetSizes(osm->lx,m,m)); 3409566063dSJacob Faibussowitsch PetscCall(VecSetType(osm->lx,vtype)); 3419566063dSJacob Faibussowitsch PetscCall(VecDuplicate(osm->lx, &osm->ly)); 3429566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(vec,osm->lis,osm->lx,isl,&osm->restriction)); 3439566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isl)); 344347574c9Seaulisa 34580ec0b7dSPatrick Sanan for (i=0; i<osm->n_local_true; ++i) { 3465b3c0d42Seaulisa ISLocalToGlobalMapping ltog; 3475b3c0d42Seaulisa IS isll; 3485b3c0d42Seaulisa const PetscInt *idx_is; 3495b3c0d42Seaulisa PetscInt *idx_lis,nout; 3505b3c0d42Seaulisa 3519566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(osm->is[i],&m)); 3529566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(osm->pmat[i],&osm->x[i],NULL)); 3539566063dSJacob Faibussowitsch PetscCall(VecDuplicate(osm->x[i],&osm->y[i])); 35455817e79SBarry Smith 355b0de9f23SBarry Smith /* generate a scatter from ly to y[i] picking all the overlapping is[i] entries */ 3569566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreateIS(osm->lis,<og)); 3579566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(osm->is[i],&m)); 3589566063dSJacob Faibussowitsch PetscCall(ISGetIndices(osm->is[i], &idx_is)); 3599566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m,&idx_lis)); 3609566063dSJacob Faibussowitsch PetscCall(ISGlobalToLocalMappingApply(ltog,IS_GTOLM_DROP,m,idx_is,&nout,idx_lis)); 3617827d75bSBarry Smith PetscCheck(nout == m,PETSC_COMM_SELF,PETSC_ERR_PLIB,"is not a subset of lis"); 3629566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(osm->is[i], &idx_is)); 3639566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF,m,idx_lis,PETSC_OWN_POINTER,&isll)); 3649566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingDestroy(<og)); 3659566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF,m,0,1,&isl)); 3669566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(osm->ly,isll,osm->y[i],isl,&osm->lrestriction[i])); 3679566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isll)); 3689566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isl)); 369910cf402Sprj- if (osm->lprolongation) { /* generate a scatter from y[i] to ly picking only the the non-overlapping is_local[i] entries */ 370d8b3b5e3Seaulisa ISLocalToGlobalMapping ltog; 371d8b3b5e3Seaulisa IS isll,isll_local; 372d8b3b5e3Seaulisa const PetscInt *idx_local; 373d8b3b5e3Seaulisa PetscInt *idx1, *idx2, nout; 374d8b3b5e3Seaulisa 3759566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(osm->is_local[i],&m_local)); 3769566063dSJacob Faibussowitsch PetscCall(ISGetIndices(osm->is_local[i], &idx_local)); 377d8b3b5e3Seaulisa 3789566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreateIS(osm->is[i],<og)); 3799566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m_local,&idx1)); 3809566063dSJacob Faibussowitsch PetscCall(ISGlobalToLocalMappingApply(ltog,IS_GTOLM_DROP,m_local,idx_local,&nout,idx1)); 3819566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingDestroy(<og)); 3827827d75bSBarry Smith PetscCheck(nout == m_local,PETSC_COMM_SELF,PETSC_ERR_PLIB,"is_local not a subset of is"); 3839566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF,m_local,idx1,PETSC_OWN_POINTER,&isll)); 384d8b3b5e3Seaulisa 3859566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingCreateIS(osm->lis,<og)); 3869566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(m_local,&idx2)); 3879566063dSJacob Faibussowitsch PetscCall(ISGlobalToLocalMappingApply(ltog,IS_GTOLM_DROP,m_local,idx_local,&nout,idx2)); 3889566063dSJacob Faibussowitsch PetscCall(ISLocalToGlobalMappingDestroy(<og)); 3897827d75bSBarry Smith PetscCheck(nout == m_local,PETSC_COMM_SELF,PETSC_ERR_PLIB,"is_local not a subset of lis"); 3909566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF,m_local,idx2,PETSC_OWN_POINTER,&isll_local)); 391d8b3b5e3Seaulisa 3929566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(osm->is_local[i], &idx_local)); 3939566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(osm->y[i],isll,osm->ly,isll_local,&osm->lprolongation[i])); 394d8b3b5e3Seaulisa 3959566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isll)); 3969566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isll_local)); 397d8b3b5e3Seaulisa } 39880ec0b7dSPatrick Sanan } 3999566063dSJacob Faibussowitsch PetscCall(VecDestroy(&vec)); 40080ec0b7dSPatrick Sanan } 40180ec0b7dSPatrick Sanan 402fb745f2cSMatthew G. Knepley if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE) { 403235cc4ceSMatthew G. Knepley IS *cis; 404235cc4ceSMatthew G. Knepley PetscInt c; 405235cc4ceSMatthew G. Knepley 4069566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(osm->n_local_true, &cis)); 407235cc4ceSMatthew G. Knepley for (c = 0; c < osm->n_local_true; ++c) cis[c] = osm->lis; 4089566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(pc->pmat, osm->n_local_true, osm->is, cis, scall, &osm->lmats)); 4099566063dSJacob Faibussowitsch PetscCall(PetscFree(cis)); 410fb745f2cSMatthew G. Knepley } 4114b9ad928SBarry Smith 4124b9ad928SBarry Smith /* Return control to the user so that the submatrices can be modified (e.g., to apply 4134b9ad928SBarry Smith different boundary conditions for the submatrices than for the global problem) */ 4149566063dSJacob Faibussowitsch PetscCall(PCModifySubMatrices(pc,osm->n_local_true,osm->is,osm->is,osm->pmat,pc->modifysubmatricesP)); 4154b9ad928SBarry Smith 41692bb6962SLisandro Dalcin /* 41792bb6962SLisandro Dalcin Loop over subdomains putting them into local ksp 41892bb6962SLisandro Dalcin */ 4199566063dSJacob Faibussowitsch PetscCall(KSPGetOptionsPrefix(osm->ksp[0],&prefix)); 42092bb6962SLisandro Dalcin for (i=0; i<osm->n_local_true; i++) { 4219566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(osm->ksp[i],osm->pmat[i],osm->pmat[i])); 4229566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(osm->pmat[i],prefix)); 42392bb6962SLisandro Dalcin if (!pc->setupcalled) { 4249566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(osm->ksp[i])); 4254b9ad928SBarry Smith } 42692bb6962SLisandro Dalcin } 4274b9ad928SBarry Smith PetscFunctionReturn(0); 4284b9ad928SBarry Smith } 4294b9ad928SBarry Smith 4306849ba73SBarry Smith static PetscErrorCode PCSetUpOnBlocks_ASM(PC pc) 4314b9ad928SBarry Smith { 4324b9ad928SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 43313f74950SBarry Smith PetscInt i; 434539c167fSBarry Smith KSPConvergedReason reason; 4354b9ad928SBarry Smith 4364b9ad928SBarry Smith PetscFunctionBegin; 4374b9ad928SBarry Smith for (i=0; i<osm->n_local_true; i++) { 4389566063dSJacob Faibussowitsch PetscCall(KSPSetUp(osm->ksp[i])); 4399566063dSJacob Faibussowitsch PetscCall(KSPGetConvergedReason(osm->ksp[i],&reason)); 440c0decd05SBarry Smith if (reason == KSP_DIVERGED_PC_FAILED) { 441261222b2SHong Zhang pc->failedreason = PC_SUBPC_ERROR; 442e0eafd54SHong Zhang } 4434b9ad928SBarry Smith } 4444b9ad928SBarry Smith PetscFunctionReturn(0); 4454b9ad928SBarry Smith } 4464b9ad928SBarry Smith 4476849ba73SBarry Smith static PetscErrorCode PCApply_ASM(PC pc,Vec x,Vec y) 4484b9ad928SBarry Smith { 4494b9ad928SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 4501dd8081eSeaulisa PetscInt i,n_local_true = osm->n_local_true; 4514b9ad928SBarry Smith ScatterMode forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE; 4524b9ad928SBarry Smith 4534b9ad928SBarry Smith PetscFunctionBegin; 4544b9ad928SBarry Smith /* 45548e38a8aSPierre Jolivet support for limiting the restriction or interpolation to only local 4564b9ad928SBarry Smith subdomain values (leaving the other values 0). 4574b9ad928SBarry Smith */ 4584b9ad928SBarry Smith if (!(osm->type & PC_ASM_RESTRICT)) { 4594b9ad928SBarry Smith forward = SCATTER_FORWARD_LOCAL; 4604b9ad928SBarry Smith /* have to zero the work RHS since scatter may leave some slots empty */ 4619566063dSJacob Faibussowitsch PetscCall(VecSet(osm->lx, 0.0)); 4624b9ad928SBarry Smith } 463347574c9Seaulisa if (!(osm->type & PC_ASM_INTERPOLATE)) { 464347574c9Seaulisa reverse = SCATTER_REVERSE_LOCAL; 465347574c9Seaulisa } 4664b9ad928SBarry Smith 467347574c9Seaulisa if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE || osm->loctype == PC_COMPOSITE_ADDITIVE) { 468b0de9f23SBarry Smith /* zero the global and the local solutions */ 4699566063dSJacob Faibussowitsch PetscCall(VecSet(y, 0.0)); 4709566063dSJacob Faibussowitsch PetscCall(VecSet(osm->ly, 0.0)); 471347574c9Seaulisa 47248e38a8aSPierre Jolivet /* copy the global RHS to local RHS including the ghost nodes */ 4739566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->restriction, x, osm->lx, INSERT_VALUES, forward)); 4749566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->restriction, x, osm->lx, INSERT_VALUES, forward)); 475347574c9Seaulisa 47648e38a8aSPierre Jolivet /* restrict local RHS to the overlapping 0-block RHS */ 4779566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward)); 4789566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward)); 479d8b3b5e3Seaulisa 48012cd4985SMatthew G. Knepley /* do the local solves */ 48112cd4985SMatthew G. Knepley for (i = 0; i < n_local_true; ++i) { 482347574c9Seaulisa 483b0de9f23SBarry Smith /* solve the overlapping i-block */ 4849566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PC_ApplyOnBlocks, osm->ksp[i], osm->x[i], osm->y[i],0)); 4859566063dSJacob Faibussowitsch PetscCall(KSPSolve(osm->ksp[i], osm->x[i], osm->y[i])); 4869566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(osm->ksp[i], pc, osm->y[i])); 4879566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PC_ApplyOnBlocks, osm->ksp[i], osm->x[i], osm->y[i], 0)); 488d8b3b5e3Seaulisa 489910cf402Sprj- if (osm->lprolongation) { /* interpolate the non-overlapping i-block solution to the local solution (only for restrictive additive) */ 4909566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward)); 4919566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward)); 49248e38a8aSPierre Jolivet } else { /* interpolate the overlapping i-block solution to the local solution */ 4939566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse)); 4949566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse)); 495d8b3b5e3Seaulisa } 496347574c9Seaulisa 497347574c9Seaulisa if (i < n_local_true-1) { 49848e38a8aSPierre Jolivet /* restrict local RHS to the overlapping (i+1)-block RHS */ 4999566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->lrestriction[i+1], osm->lx, osm->x[i+1], INSERT_VALUES, forward)); 5009566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->lrestriction[i+1], osm->lx, osm->x[i+1], INSERT_VALUES, forward)); 501347574c9Seaulisa 502347574c9Seaulisa if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE) { 5038966356dSPierre Jolivet /* update the overlapping (i+1)-block RHS using the current local solution */ 5049566063dSJacob Faibussowitsch PetscCall(MatMult(osm->lmats[i+1], osm->ly, osm->y[i+1])); 5059566063dSJacob Faibussowitsch PetscCall(VecAXPBY(osm->x[i+1],-1.,1., osm->y[i+1])); 5067c3d802fSMatthew G. Knepley } 50712cd4985SMatthew G. Knepley } 50812cd4985SMatthew G. Knepley } 50948e38a8aSPierre Jolivet /* add the local solution to the global solution including the ghost nodes */ 5109566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->restriction, osm->ly, y, ADD_VALUES, reverse)); 5119566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->restriction, osm->ly, y, ADD_VALUES, reverse)); 51298921bdaSJacob Faibussowitsch } else SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Invalid local composition type: %s", PCCompositeTypes[osm->loctype]); 5134b9ad928SBarry Smith PetscFunctionReturn(0); 5144b9ad928SBarry Smith } 5154b9ad928SBarry Smith 51648e38a8aSPierre Jolivet static PetscErrorCode PCMatApply_ASM(PC pc,Mat X,Mat Y) 51748e38a8aSPierre Jolivet { 51848e38a8aSPierre Jolivet PC_ASM *osm = (PC_ASM*)pc->data; 51948e38a8aSPierre Jolivet Mat Z,W; 52048e38a8aSPierre Jolivet Vec x; 52148e38a8aSPierre Jolivet PetscInt i,m,N; 52248e38a8aSPierre Jolivet ScatterMode forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE; 52348e38a8aSPierre Jolivet 52448e38a8aSPierre Jolivet PetscFunctionBegin; 5257827d75bSBarry Smith PetscCheck(osm->n_local_true <= 1,PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Not yet implemented"); 52648e38a8aSPierre Jolivet /* 52748e38a8aSPierre Jolivet support for limiting the restriction or interpolation to only local 52848e38a8aSPierre Jolivet subdomain values (leaving the other values 0). 52948e38a8aSPierre Jolivet */ 53048e38a8aSPierre Jolivet if (!(osm->type & PC_ASM_RESTRICT)) { 53148e38a8aSPierre Jolivet forward = SCATTER_FORWARD_LOCAL; 53248e38a8aSPierre Jolivet /* have to zero the work RHS since scatter may leave some slots empty */ 5339566063dSJacob Faibussowitsch PetscCall(VecSet(osm->lx, 0.0)); 53448e38a8aSPierre Jolivet } 53548e38a8aSPierre Jolivet if (!(osm->type & PC_ASM_INTERPOLATE)) { 53648e38a8aSPierre Jolivet reverse = SCATTER_REVERSE_LOCAL; 53748e38a8aSPierre Jolivet } 5389566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(osm->x[0], &m)); 5399566063dSJacob Faibussowitsch PetscCall(MatGetSize(X, NULL, &N)); 5409566063dSJacob Faibussowitsch PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, m, N, NULL, &Z)); 54148e38a8aSPierre Jolivet if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE || osm->loctype == PC_COMPOSITE_ADDITIVE) { 54248e38a8aSPierre Jolivet /* zero the global and the local solutions */ 5439566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(Y)); 5449566063dSJacob Faibussowitsch PetscCall(VecSet(osm->ly, 0.0)); 54548e38a8aSPierre Jolivet 54648e38a8aSPierre Jolivet for (i = 0; i < N; ++i) { 5479566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumnVecRead(X, i, &x)); 54848e38a8aSPierre Jolivet /* copy the global RHS to local RHS including the ghost nodes */ 5499566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->restriction, x, osm->lx, INSERT_VALUES, forward)); 5509566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->restriction, x, osm->lx, INSERT_VALUES, forward)); 5519566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumnVecRead(X, i, &x)); 55248e38a8aSPierre Jolivet 5539566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumnVecWrite(Z, i, &x)); 55448e38a8aSPierre Jolivet /* restrict local RHS to the overlapping 0-block RHS */ 5559566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->lrestriction[0], osm->lx, x, INSERT_VALUES, forward)); 5569566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->lrestriction[0], osm->lx, x, INSERT_VALUES, forward)); 5579566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumnVecWrite(Z, i, &x)); 55848e38a8aSPierre Jolivet } 5599566063dSJacob Faibussowitsch PetscCall(MatCreateSeqDense(PETSC_COMM_SELF, m, N, NULL, &W)); 56048e38a8aSPierre Jolivet /* solve the overlapping 0-block */ 5619566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PC_ApplyOnBlocks, osm->ksp[0], Z, W, 0)); 5629566063dSJacob Faibussowitsch PetscCall(KSPMatSolve(osm->ksp[0], Z, W)); 5639566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(osm->ksp[0], pc, NULL)); 5649566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PC_ApplyOnBlocks, osm->ksp[0], Z, W,0)); 5659566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Z)); 56648e38a8aSPierre Jolivet 56748e38a8aSPierre Jolivet for (i = 0; i < N; ++i) { 5689566063dSJacob Faibussowitsch PetscCall(VecSet(osm->ly, 0.0)); 5699566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumnVecRead(W, i, &x)); 57048e38a8aSPierre Jolivet if (osm->lprolongation) { /* interpolate the non-overlapping 0-block solution to the local solution (only for restrictive additive) */ 5719566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->lprolongation[0], x, osm->ly, ADD_VALUES, forward)); 5729566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->lprolongation[0], x, osm->ly, ADD_VALUES, forward)); 57348e38a8aSPierre Jolivet } else { /* interpolate the overlapping 0-block solution to the local solution */ 5749566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->lrestriction[0], x, osm->ly, ADD_VALUES, reverse)); 5759566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->lrestriction[0], x, osm->ly, ADD_VALUES, reverse)); 57648e38a8aSPierre Jolivet } 5779566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumnVecRead(W, i, &x)); 57848e38a8aSPierre Jolivet 5799566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumnVecWrite(Y, i, &x)); 58048e38a8aSPierre Jolivet /* add the local solution to the global solution including the ghost nodes */ 5819566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->restriction, osm->ly, x, ADD_VALUES, reverse)); 5829566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->restriction, osm->ly, x, ADD_VALUES, reverse)); 5839566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumnVecWrite(Y, i, &x)); 58448e38a8aSPierre Jolivet } 5859566063dSJacob Faibussowitsch PetscCall(MatDestroy(&W)); 58698921bdaSJacob Faibussowitsch } else SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Invalid local composition type: %s", PCCompositeTypes[osm->loctype]); 58748e38a8aSPierre Jolivet PetscFunctionReturn(0); 58848e38a8aSPierre Jolivet } 58948e38a8aSPierre Jolivet 5906849ba73SBarry Smith static PetscErrorCode PCApplyTranspose_ASM(PC pc,Vec x,Vec y) 5914b9ad928SBarry Smith { 5924b9ad928SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 5931dd8081eSeaulisa PetscInt i,n_local_true = osm->n_local_true; 5944b9ad928SBarry Smith ScatterMode forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE; 5954b9ad928SBarry Smith 5964b9ad928SBarry Smith PetscFunctionBegin; 5974b9ad928SBarry Smith /* 5984b9ad928SBarry Smith Support for limiting the restriction or interpolation to only local 5994b9ad928SBarry Smith subdomain values (leaving the other values 0). 6004b9ad928SBarry Smith 6014b9ad928SBarry Smith Note: these are reversed from the PCApply_ASM() because we are applying the 6024b9ad928SBarry Smith transpose of the three terms 6034b9ad928SBarry Smith */ 604d8b3b5e3Seaulisa 6054b9ad928SBarry Smith if (!(osm->type & PC_ASM_INTERPOLATE)) { 6064b9ad928SBarry Smith forward = SCATTER_FORWARD_LOCAL; 6074b9ad928SBarry Smith /* have to zero the work RHS since scatter may leave some slots empty */ 6089566063dSJacob Faibussowitsch PetscCall(VecSet(osm->lx, 0.0)); 6094b9ad928SBarry Smith } 6102fa5cd67SKarl Rupp if (!(osm->type & PC_ASM_RESTRICT)) reverse = SCATTER_REVERSE_LOCAL; 6114b9ad928SBarry Smith 612b0de9f23SBarry Smith /* zero the global and the local solutions */ 6139566063dSJacob Faibussowitsch PetscCall(VecSet(y, 0.0)); 6149566063dSJacob Faibussowitsch PetscCall(VecSet(osm->ly, 0.0)); 615d8b3b5e3Seaulisa 616b0de9f23SBarry Smith /* Copy the global RHS to local RHS including the ghost nodes */ 6179566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->restriction, x, osm->lx, INSERT_VALUES, forward)); 6189566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->restriction, x, osm->lx, INSERT_VALUES, forward)); 619d8b3b5e3Seaulisa 620b0de9f23SBarry Smith /* Restrict local RHS to the overlapping 0-block RHS */ 6219566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward)); 6229566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward)); 623d8b3b5e3Seaulisa 6244b9ad928SBarry Smith /* do the local solves */ 625d8b3b5e3Seaulisa for (i = 0; i < n_local_true; ++i) { 626d8b3b5e3Seaulisa 627b0de9f23SBarry Smith /* solve the overlapping i-block */ 6289566063dSJacob Faibussowitsch PetscCall(PetscLogEventBegin(PC_ApplyOnBlocks,osm->ksp[i],osm->x[i],osm->y[i],0)); 6299566063dSJacob Faibussowitsch PetscCall(KSPSolveTranspose(osm->ksp[i], osm->x[i], osm->y[i])); 6309566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(osm->ksp[i],pc,osm->y[i])); 6319566063dSJacob Faibussowitsch PetscCall(PetscLogEventEnd(PC_ApplyOnBlocks,osm->ksp[i],osm->x[i],osm->y[i],0)); 632d8b3b5e3Seaulisa 633910cf402Sprj- if (osm->lprolongation) { /* interpolate the non-overlapping i-block solution to the local solution */ 6349566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward)); 6359566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward)); 636910cf402Sprj- } else { /* interpolate the overlapping i-block solution to the local solution */ 6379566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse)); 6389566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse)); 6394b9ad928SBarry Smith } 640d8b3b5e3Seaulisa 641d8b3b5e3Seaulisa if (i < n_local_true-1) { 642b0de9f23SBarry Smith /* Restrict local RHS to the overlapping (i+1)-block RHS */ 6439566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->lrestriction[i+1], osm->lx, osm->x[i+1], INSERT_VALUES, forward)); 6449566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->lrestriction[i+1], osm->lx, osm->x[i+1], INSERT_VALUES, forward)); 6454b9ad928SBarry Smith } 6464b9ad928SBarry Smith } 647b0de9f23SBarry Smith /* Add the local solution to the global solution including the ghost nodes */ 6489566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->restriction, osm->ly, y, ADD_VALUES, reverse)); 6499566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->restriction, osm->ly, y, ADD_VALUES, reverse)); 650d8b3b5e3Seaulisa 6514b9ad928SBarry Smith PetscFunctionReturn(0); 652d8b3b5e3Seaulisa 6534b9ad928SBarry Smith } 6544b9ad928SBarry Smith 655e91c6855SBarry Smith static PetscErrorCode PCReset_ASM(PC pc) 6564b9ad928SBarry Smith { 6574b9ad928SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 65813f74950SBarry Smith PetscInt i; 6594b9ad928SBarry Smith 6604b9ad928SBarry Smith PetscFunctionBegin; 66192bb6962SLisandro Dalcin if (osm->ksp) { 66292bb6962SLisandro Dalcin for (i=0; i<osm->n_local_true; i++) { 6639566063dSJacob Faibussowitsch PetscCall(KSPReset(osm->ksp[i])); 66492bb6962SLisandro Dalcin } 66592bb6962SLisandro Dalcin } 666e09e08ccSBarry Smith if (osm->pmat) { 66792bb6962SLisandro Dalcin if (osm->n_local_true > 0) { 6689566063dSJacob Faibussowitsch PetscCall(MatDestroySubMatrices(osm->n_local_true,&osm->pmat)); 66992bb6962SLisandro Dalcin } 67092bb6962SLisandro Dalcin } 6711dd8081eSeaulisa if (osm->lrestriction) { 6729566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&osm->restriction)); 6731dd8081eSeaulisa for (i=0; i<osm->n_local_true; i++) { 6749566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&osm->lrestriction[i])); 6759566063dSJacob Faibussowitsch if (osm->lprolongation) PetscCall(VecScatterDestroy(&osm->lprolongation[i])); 6769566063dSJacob Faibussowitsch PetscCall(VecDestroy(&osm->x[i])); 6779566063dSJacob Faibussowitsch PetscCall(VecDestroy(&osm->y[i])); 6784b9ad928SBarry Smith } 6799566063dSJacob Faibussowitsch PetscCall(PetscFree(osm->lrestriction)); 6809566063dSJacob Faibussowitsch if (osm->lprolongation) PetscCall(PetscFree(osm->lprolongation)); 6819566063dSJacob Faibussowitsch PetscCall(PetscFree(osm->x)); 6829566063dSJacob Faibussowitsch PetscCall(PetscFree(osm->y)); 6831dd8081eSeaulisa 68492bb6962SLisandro Dalcin } 6859566063dSJacob Faibussowitsch PetscCall(PCASMDestroySubdomains(osm->n_local_true,osm->is,osm->is_local)); 6869566063dSJacob Faibussowitsch PetscCall(ISDestroy(&osm->lis)); 6879566063dSJacob Faibussowitsch PetscCall(VecDestroy(&osm->lx)); 6889566063dSJacob Faibussowitsch PetscCall(VecDestroy(&osm->ly)); 689347574c9Seaulisa if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE) { 6909566063dSJacob Faibussowitsch PetscCall(MatDestroyMatrices(osm->n_local_true, &osm->lmats)); 691fb745f2cSMatthew G. Knepley } 6922fa5cd67SKarl Rupp 6939566063dSJacob Faibussowitsch PetscCall(PetscFree(osm->sub_mat_type)); 69480ec0b7dSPatrick Sanan 6950a545947SLisandro Dalcin osm->is = NULL; 6960a545947SLisandro Dalcin osm->is_local = NULL; 697e91c6855SBarry Smith PetscFunctionReturn(0); 698e91c6855SBarry Smith } 699e91c6855SBarry Smith 700e91c6855SBarry Smith static PetscErrorCode PCDestroy_ASM(PC pc) 701e91c6855SBarry Smith { 702e91c6855SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 703e91c6855SBarry Smith PetscInt i; 704e91c6855SBarry Smith 705e91c6855SBarry Smith PetscFunctionBegin; 7069566063dSJacob Faibussowitsch PetscCall(PCReset_ASM(pc)); 707e91c6855SBarry Smith if (osm->ksp) { 708e91c6855SBarry Smith for (i=0; i<osm->n_local_true; i++) { 7099566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&osm->ksp[i])); 710e91c6855SBarry Smith } 7119566063dSJacob Faibussowitsch PetscCall(PetscFree(osm->ksp)); 712e91c6855SBarry Smith } 7139566063dSJacob Faibussowitsch PetscCall(PetscFree(pc->data)); 71496322394SPierre Jolivet 7159566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCASMSetLocalSubdomains_C",NULL)); 7169566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCASMSetTotalSubdomains_C",NULL)); 7179566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCASMSetOverlap_C",NULL)); 7189566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCASMSetType_C",NULL)); 7199566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCASMGetType_C",NULL)); 7209566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCASMSetLocalType_C",NULL)); 7219566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCASMGetLocalType_C",NULL)); 7229566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCASMSetSortIndices_C",NULL)); 7239566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCASMGetSubKSP_C",NULL)); 7249566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCASMGetSubMatType_C",NULL)); 7259566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCASMSetSubMatType_C",NULL)); 7264b9ad928SBarry Smith PetscFunctionReturn(0); 7274b9ad928SBarry Smith } 7284b9ad928SBarry Smith 7294416b707SBarry Smith static PetscErrorCode PCSetFromOptions_ASM(PetscOptionItems *PetscOptionsObject,PC pc) 7304b9ad928SBarry Smith { 7314b9ad928SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 7329dcbbd2bSBarry Smith PetscInt blocks,ovl; 73357501b6eSBarry Smith PetscBool flg; 73492bb6962SLisandro Dalcin PCASMType asmtype; 73512cd4985SMatthew G. Knepley PCCompositeType loctype; 73680ec0b7dSPatrick Sanan char sub_mat_type[256]; 7374b9ad928SBarry Smith 7384b9ad928SBarry Smith PetscFunctionBegin; 739d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject,"Additive Schwarz options"); 7409566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_asm_dm_subdomains","Use DMCreateDomainDecomposition() to define subdomains","PCASMSetDMSubdomains",osm->dm_subdomains,&osm->dm_subdomains,&flg)); 7419566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_asm_blocks","Number of subdomains","PCASMSetTotalSubdomains",osm->n,&blocks,&flg)); 74265db9045SDmitry Karpeev if (flg) { 7439566063dSJacob Faibussowitsch PetscCall(PCASMSetTotalSubdomains(pc,blocks,NULL,NULL)); 744d709ea83SDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 74565db9045SDmitry Karpeev } 7469566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_asm_local_blocks","Number of local subdomains","PCASMSetLocalSubdomains",osm->n_local_true,&blocks,&flg)); 747342c94f9SMatthew G. Knepley if (flg) { 7489566063dSJacob Faibussowitsch PetscCall(PCASMSetLocalSubdomains(pc,blocks,NULL,NULL)); 749342c94f9SMatthew G. Knepley osm->dm_subdomains = PETSC_FALSE; 750342c94f9SMatthew G. Knepley } 7519566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_asm_overlap","Number of grid points overlap","PCASMSetOverlap",osm->overlap,&ovl,&flg)); 75265db9045SDmitry Karpeev if (flg) { 7539566063dSJacob Faibussowitsch PetscCall(PCASMSetOverlap(pc,ovl)); 754d709ea83SDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 75565db9045SDmitry Karpeev } 75690d69ab7SBarry Smith flg = PETSC_FALSE; 7579566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-pc_asm_type","Type of restriction/extension","PCASMSetType",PCASMTypes,(PetscEnum)osm->type,(PetscEnum*)&asmtype,&flg)); 7589566063dSJacob Faibussowitsch if (flg) PetscCall(PCASMSetType(pc,asmtype)); 75912cd4985SMatthew G. Knepley flg = PETSC_FALSE; 7609566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-pc_asm_local_type","Type of local solver composition","PCASMSetLocalType",PCCompositeTypes,(PetscEnum)osm->loctype,(PetscEnum*)&loctype,&flg)); 7619566063dSJacob Faibussowitsch if (flg) PetscCall(PCASMSetLocalType(pc,loctype)); 7629566063dSJacob Faibussowitsch PetscCall(PetscOptionsFList("-pc_asm_sub_mat_type","Subsolve Matrix Type","PCASMSetSubMatType",MatList,NULL,sub_mat_type,256,&flg)); 76380ec0b7dSPatrick Sanan if (flg) { 7649566063dSJacob Faibussowitsch PetscCall(PCASMSetSubMatType(pc,sub_mat_type)); 76580ec0b7dSPatrick Sanan } 766d0609cedSBarry Smith PetscOptionsHeadEnd(); 7674b9ad928SBarry Smith PetscFunctionReturn(0); 7684b9ad928SBarry Smith } 7694b9ad928SBarry Smith 7704b9ad928SBarry Smith /*------------------------------------------------------------------------------------*/ 7714b9ad928SBarry Smith 7721e6b0712SBarry Smith static PetscErrorCode PCASMSetLocalSubdomains_ASM(PC pc,PetscInt n,IS is[],IS is_local[]) 7734b9ad928SBarry Smith { 7744b9ad928SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 77592bb6962SLisandro Dalcin PetscInt i; 7764b9ad928SBarry Smith 7774b9ad928SBarry Smith PetscFunctionBegin; 77863a3b9bcSJacob Faibussowitsch PetscCheck(n >= 1,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Each process must have 1 or more blocks, n = %" PetscInt_FMT,n); 7797827d75bSBarry Smith PetscCheck(!pc->setupcalled || (n == osm->n_local_true && !is),PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCASMSetLocalSubdomains() should be called before calling PCSetUp()."); 780e7e72b3dSBarry Smith 7814b9ad928SBarry Smith if (!pc->setupcalled) { 78292bb6962SLisandro Dalcin if (is) { 7839566063dSJacob Faibussowitsch for (i=0; i<n; i++) PetscCall(PetscObjectReference((PetscObject)is[i])); 78492bb6962SLisandro Dalcin } 785832fc9a5SMatthew Knepley if (is_local) { 7869566063dSJacob Faibussowitsch for (i=0; i<n; i++) PetscCall(PetscObjectReference((PetscObject)is_local[i])); 787832fc9a5SMatthew Knepley } 7889566063dSJacob Faibussowitsch PetscCall(PCASMDestroySubdomains(osm->n_local_true,osm->is,osm->is_local)); 7892fa5cd67SKarl Rupp 7904b9ad928SBarry Smith osm->n_local_true = n; 7910a545947SLisandro Dalcin osm->is = NULL; 7920a545947SLisandro Dalcin osm->is_local = NULL; 79392bb6962SLisandro Dalcin if (is) { 7949566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n,&osm->is)); 7952fa5cd67SKarl Rupp for (i=0; i<n; i++) osm->is[i] = is[i]; 7963d03bb48SJed Brown /* Flag indicating that the user has set overlapping subdomains so PCASM should not increase their size. */ 7973d03bb48SJed Brown osm->overlap = -1; 79892bb6962SLisandro Dalcin } 7992b691e39SMatthew Knepley if (is_local) { 8009566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n,&osm->is_local)); 8012fa5cd67SKarl Rupp for (i=0; i<n; i++) osm->is_local[i] = is_local[i]; 802a35b7b57SDmitry Karpeev if (!is) { 8039566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(osm->n_local_true,&osm->is)); 804a35b7b57SDmitry Karpeev for (i=0; i<osm->n_local_true; i++) { 805a35b7b57SDmitry Karpeev if (osm->overlap > 0) { /* With positive overlap, osm->is[i] will be modified */ 8069566063dSJacob Faibussowitsch PetscCall(ISDuplicate(osm->is_local[i],&osm->is[i])); 8079566063dSJacob Faibussowitsch PetscCall(ISCopy(osm->is_local[i],osm->is[i])); 808a35b7b57SDmitry Karpeev } else { 8099566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)osm->is_local[i])); 810a35b7b57SDmitry Karpeev osm->is[i] = osm->is_local[i]; 811a35b7b57SDmitry Karpeev } 812a35b7b57SDmitry Karpeev } 813a35b7b57SDmitry Karpeev } 8142b691e39SMatthew Knepley } 8154b9ad928SBarry Smith } 8164b9ad928SBarry Smith PetscFunctionReturn(0); 8174b9ad928SBarry Smith } 8184b9ad928SBarry Smith 8191e6b0712SBarry Smith static PetscErrorCode PCASMSetTotalSubdomains_ASM(PC pc,PetscInt N,IS *is,IS *is_local) 8204b9ad928SBarry Smith { 8214b9ad928SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 82213f74950SBarry Smith PetscMPIInt rank,size; 82378904715SLisandro Dalcin PetscInt n; 8244b9ad928SBarry Smith 8254b9ad928SBarry Smith PetscFunctionBegin; 82663a3b9bcSJacob Faibussowitsch PetscCheck(N >= 1,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Number of total blocks must be > 0, N = %" PetscInt_FMT,N); 8277827d75bSBarry 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."); 8284b9ad928SBarry Smith 8294b9ad928SBarry Smith /* 830880770e9SJed Brown Split the subdomains equally among all processors 8314b9ad928SBarry Smith */ 8329566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank)); 8339566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size)); 8344b9ad928SBarry Smith n = N/size + ((N % size) > rank); 83563a3b9bcSJacob 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); 8367827d75bSBarry Smith PetscCheck(!pc->setupcalled || n == osm->n_local_true,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"PCASMSetTotalSubdomains() should be called before PCSetUp()."); 8374b9ad928SBarry Smith if (!pc->setupcalled) { 8389566063dSJacob Faibussowitsch PetscCall(PCASMDestroySubdomains(osm->n_local_true,osm->is,osm->is_local)); 8392fa5cd67SKarl Rupp 8404b9ad928SBarry Smith osm->n_local_true = n; 8410a545947SLisandro Dalcin osm->is = NULL; 8420a545947SLisandro Dalcin osm->is_local = NULL; 8434b9ad928SBarry Smith } 8444b9ad928SBarry Smith PetscFunctionReturn(0); 8454b9ad928SBarry Smith } 8464b9ad928SBarry Smith 8471e6b0712SBarry Smith static PetscErrorCode PCASMSetOverlap_ASM(PC pc,PetscInt ovl) 8484b9ad928SBarry Smith { 84992bb6962SLisandro Dalcin PC_ASM *osm = (PC_ASM*)pc->data; 8504b9ad928SBarry Smith 8514b9ad928SBarry Smith PetscFunctionBegin; 8527827d75bSBarry Smith PetscCheck(ovl >= 0,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap value requested"); 8537827d75bSBarry Smith PetscCheck(!pc->setupcalled || ovl == osm->overlap,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCASMSetOverlap() should be called before PCSetUp()."); 8542fa5cd67SKarl Rupp if (!pc->setupcalled) osm->overlap = ovl; 8554b9ad928SBarry Smith PetscFunctionReturn(0); 8564b9ad928SBarry Smith } 8574b9ad928SBarry Smith 8581e6b0712SBarry Smith static PetscErrorCode PCASMSetType_ASM(PC pc,PCASMType type) 8594b9ad928SBarry Smith { 86092bb6962SLisandro Dalcin PC_ASM *osm = (PC_ASM*)pc->data; 8614b9ad928SBarry Smith 8624b9ad928SBarry Smith PetscFunctionBegin; 8634b9ad928SBarry Smith osm->type = type; 864bf108f30SBarry Smith osm->type_set = PETSC_TRUE; 8654b9ad928SBarry Smith PetscFunctionReturn(0); 8664b9ad928SBarry Smith } 8674b9ad928SBarry Smith 868c60c7ad4SBarry Smith static PetscErrorCode PCASMGetType_ASM(PC pc,PCASMType *type) 869c60c7ad4SBarry Smith { 870c60c7ad4SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 871c60c7ad4SBarry Smith 872c60c7ad4SBarry Smith PetscFunctionBegin; 873c60c7ad4SBarry Smith *type = osm->type; 874c60c7ad4SBarry Smith PetscFunctionReturn(0); 875c60c7ad4SBarry Smith } 876c60c7ad4SBarry Smith 87712cd4985SMatthew G. Knepley static PetscErrorCode PCASMSetLocalType_ASM(PC pc, PCCompositeType type) 87812cd4985SMatthew G. Knepley { 87912cd4985SMatthew G. Knepley PC_ASM *osm = (PC_ASM *) pc->data; 88012cd4985SMatthew G. Knepley 88112cd4985SMatthew G. Knepley PetscFunctionBegin; 8827827d75bSBarry 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"); 88312cd4985SMatthew G. Knepley osm->loctype = type; 88412cd4985SMatthew G. Knepley PetscFunctionReturn(0); 88512cd4985SMatthew G. Knepley } 88612cd4985SMatthew G. Knepley 88712cd4985SMatthew G. Knepley static PetscErrorCode PCASMGetLocalType_ASM(PC pc, PCCompositeType *type) 88812cd4985SMatthew G. Knepley { 88912cd4985SMatthew G. Knepley PC_ASM *osm = (PC_ASM *) pc->data; 89012cd4985SMatthew G. Knepley 89112cd4985SMatthew G. Knepley PetscFunctionBegin; 89212cd4985SMatthew G. Knepley *type = osm->loctype; 89312cd4985SMatthew G. Knepley PetscFunctionReturn(0); 89412cd4985SMatthew G. Knepley } 89512cd4985SMatthew G. Knepley 8961e6b0712SBarry Smith static PetscErrorCode PCASMSetSortIndices_ASM(PC pc,PetscBool doSort) 8976ed231c7SMatthew Knepley { 8986ed231c7SMatthew Knepley PC_ASM *osm = (PC_ASM*)pc->data; 8996ed231c7SMatthew Knepley 9006ed231c7SMatthew Knepley PetscFunctionBegin; 9016ed231c7SMatthew Knepley osm->sort_indices = doSort; 9026ed231c7SMatthew Knepley PetscFunctionReturn(0); 9036ed231c7SMatthew Knepley } 9046ed231c7SMatthew Knepley 9051e6b0712SBarry Smith static PetscErrorCode PCASMGetSubKSP_ASM(PC pc,PetscInt *n_local,PetscInt *first_local,KSP **ksp) 9064b9ad928SBarry Smith { 90792bb6962SLisandro Dalcin PC_ASM *osm = (PC_ASM*)pc->data; 9084b9ad928SBarry Smith 9094b9ad928SBarry Smith PetscFunctionBegin; 9107827d75bSBarry 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"); 9114b9ad928SBarry Smith 9122fa5cd67SKarl Rupp if (n_local) *n_local = osm->n_local_true; 91392bb6962SLisandro Dalcin if (first_local) { 9149566063dSJacob Faibussowitsch PetscCallMPI(MPI_Scan(&osm->n_local_true,first_local,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)pc))); 91592bb6962SLisandro Dalcin *first_local -= osm->n_local_true; 91692bb6962SLisandro Dalcin } 917ed155784SPierre Jolivet if (ksp) *ksp = osm->ksp; 9184b9ad928SBarry Smith PetscFunctionReturn(0); 9194b9ad928SBarry Smith } 9204b9ad928SBarry Smith 92180ec0b7dSPatrick Sanan static PetscErrorCode PCASMGetSubMatType_ASM(PC pc,MatType *sub_mat_type) 92280ec0b7dSPatrick Sanan { 92380ec0b7dSPatrick Sanan PC_ASM *osm = (PC_ASM*)pc->data; 92480ec0b7dSPatrick Sanan 92580ec0b7dSPatrick Sanan PetscFunctionBegin; 92680ec0b7dSPatrick Sanan PetscValidHeaderSpecific(pc,PC_CLASSID,1); 92780ec0b7dSPatrick Sanan PetscValidPointer(sub_mat_type,2); 92880ec0b7dSPatrick Sanan *sub_mat_type = osm->sub_mat_type; 92980ec0b7dSPatrick Sanan PetscFunctionReturn(0); 93080ec0b7dSPatrick Sanan } 93180ec0b7dSPatrick Sanan 93280ec0b7dSPatrick Sanan static PetscErrorCode PCASMSetSubMatType_ASM(PC pc,MatType sub_mat_type) 93380ec0b7dSPatrick Sanan { 93480ec0b7dSPatrick Sanan PC_ASM *osm = (PC_ASM*)pc->data; 93580ec0b7dSPatrick Sanan 93680ec0b7dSPatrick Sanan PetscFunctionBegin; 93780ec0b7dSPatrick Sanan PetscValidHeaderSpecific(pc,PC_CLASSID,1); 9389566063dSJacob Faibussowitsch PetscCall(PetscFree(osm->sub_mat_type)); 9399566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(sub_mat_type,(char**)&osm->sub_mat_type)); 94080ec0b7dSPatrick Sanan PetscFunctionReturn(0); 94180ec0b7dSPatrick Sanan } 94280ec0b7dSPatrick Sanan 9434b9ad928SBarry Smith /*@C 9441093a601SBarry Smith PCASMSetLocalSubdomains - Sets the local subdomains (for this processor only) for the additive Schwarz preconditioner. 9454b9ad928SBarry Smith 946d083f849SBarry Smith Collective on pc 9474b9ad928SBarry Smith 9484b9ad928SBarry Smith Input Parameters: 9494b9ad928SBarry Smith + pc - the preconditioner context 9504b9ad928SBarry Smith . n - the number of subdomains for this processor (default value = 1) 9518c03b21aSDmitry Karpeev . is - the index set that defines the subdomains for this processor 9520298fd71SBarry Smith (or NULL for PETSc to determine subdomains) 953f1ee410cSBarry 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 954f1ee410cSBarry Smith (or NULL to not provide these) 9554b9ad928SBarry Smith 956342c94f9SMatthew G. Knepley Options Database Key: 957342c94f9SMatthew G. Knepley To set the total number of subdomain blocks rather than specify the 958342c94f9SMatthew G. Knepley index sets, use the option 959342c94f9SMatthew G. Knepley . -pc_asm_local_blocks <blks> - Sets local blocks 960342c94f9SMatthew G. Knepley 9614b9ad928SBarry Smith Notes: 9621093a601SBarry Smith The IS numbering is in the parallel, global numbering of the vector for both is and is_local 9634b9ad928SBarry Smith 9644b9ad928SBarry Smith By default the ASM preconditioner uses 1 block per processor. 9654b9ad928SBarry Smith 9664b9ad928SBarry Smith Use PCASMSetTotalSubdomains() to set the subdomains for all processors. 9674b9ad928SBarry Smith 968f1ee410cSBarry Smith If is_local is provided and PCASMType is PC_ASM_RESTRICT then the solution only over the is_local region is interpolated 969f1ee410cSBarry Smith back to form the global solution (this is the standard restricted additive Schwarz method) 970f1ee410cSBarry Smith 971f1ee410cSBarry 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 972f1ee410cSBarry Smith no code to handle that case. 973f1ee410cSBarry Smith 9744b9ad928SBarry Smith Level: advanced 9754b9ad928SBarry Smith 976db781477SPatrick Sanan .seealso: `PCASMSetTotalSubdomains()`, `PCASMSetOverlap()`, `PCASMGetSubKSP()`, 977db781477SPatrick Sanan `PCASMCreateSubdomains2D()`, `PCASMGetLocalSubdomains()`, `PCASMType`, `PCASMSetType()` 9784b9ad928SBarry Smith @*/ 9797087cfbeSBarry Smith PetscErrorCode PCASMSetLocalSubdomains(PC pc,PetscInt n,IS is[],IS is_local[]) 9804b9ad928SBarry Smith { 9814b9ad928SBarry Smith PetscFunctionBegin; 9820700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 983cac4c232SBarry Smith PetscTryMethod(pc,"PCASMSetLocalSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,n,is,is_local)); 9844b9ad928SBarry Smith PetscFunctionReturn(0); 9854b9ad928SBarry Smith } 9864b9ad928SBarry Smith 9874b9ad928SBarry Smith /*@C 988feb221f8SDmitry Karpeev PCASMSetTotalSubdomains - Sets the subdomains for all processors for the 9894b9ad928SBarry Smith additive Schwarz preconditioner. Either all or no processors in the 9904b9ad928SBarry Smith PC communicator must call this routine, with the same index sets. 9914b9ad928SBarry Smith 992d083f849SBarry Smith Collective on pc 9934b9ad928SBarry Smith 9944b9ad928SBarry Smith Input Parameters: 9954b9ad928SBarry Smith + pc - the preconditioner context 996feb221f8SDmitry Karpeev . N - the number of subdomains for all processors 997feb221f8SDmitry Karpeev . is - the index sets that define the subdomains for all processors 998dfaa0c5fSBarry Smith (or NULL to ask PETSc to determine the subdomains) 9992b691e39SMatthew Knepley - is_local - the index sets that define the local part of the subdomains for this processor 1000f1ee410cSBarry Smith (or NULL to not provide this information) 10014b9ad928SBarry Smith 10024b9ad928SBarry Smith Options Database Key: 10034b9ad928SBarry Smith To set the total number of subdomain blocks rather than specify the 10044b9ad928SBarry Smith index sets, use the option 10054b9ad928SBarry Smith . -pc_asm_blocks <blks> - Sets total blocks 10064b9ad928SBarry Smith 10074b9ad928SBarry Smith Notes: 1008f1ee410cSBarry Smith Currently you cannot use this to set the actual subdomains with the argument is or is_local. 10094b9ad928SBarry Smith 10104b9ad928SBarry Smith By default the ASM preconditioner uses 1 block per processor. 10114b9ad928SBarry Smith 10124b9ad928SBarry Smith These index sets cannot be destroyed until after completion of the 10134b9ad928SBarry Smith linear solves for which the ASM preconditioner is being used. 10144b9ad928SBarry Smith 10154b9ad928SBarry Smith Use PCASMSetLocalSubdomains() to set local subdomains. 10164b9ad928SBarry Smith 10171093a601SBarry Smith The IS numbering is in the parallel, global numbering of the vector for both is and is_local 10181093a601SBarry Smith 10194b9ad928SBarry Smith Level: advanced 10204b9ad928SBarry Smith 1021db781477SPatrick Sanan .seealso: `PCASMSetLocalSubdomains()`, `PCASMSetOverlap()`, `PCASMGetSubKSP()`, 1022db781477SPatrick Sanan `PCASMCreateSubdomains2D()` 10234b9ad928SBarry Smith @*/ 10247087cfbeSBarry Smith PetscErrorCode PCASMSetTotalSubdomains(PC pc,PetscInt N,IS is[],IS is_local[]) 10254b9ad928SBarry Smith { 10264b9ad928SBarry Smith PetscFunctionBegin; 10270700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1028cac4c232SBarry Smith PetscTryMethod(pc,"PCASMSetTotalSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,N,is,is_local)); 10294b9ad928SBarry Smith PetscFunctionReturn(0); 10304b9ad928SBarry Smith } 10314b9ad928SBarry Smith 10324b9ad928SBarry Smith /*@ 10334b9ad928SBarry Smith PCASMSetOverlap - Sets the overlap between a pair of subdomains for the 10344b9ad928SBarry Smith additive Schwarz preconditioner. Either all or no processors in the 1035f1ee410cSBarry Smith PC communicator must call this routine. 10364b9ad928SBarry Smith 1037d083f849SBarry Smith Logically Collective on pc 10384b9ad928SBarry Smith 10394b9ad928SBarry Smith Input Parameters: 10404b9ad928SBarry Smith + pc - the preconditioner context 10414b9ad928SBarry Smith - ovl - the amount of overlap between subdomains (ovl >= 0, default value = 1) 10424b9ad928SBarry Smith 10434b9ad928SBarry Smith Options Database Key: 10444b9ad928SBarry Smith . -pc_asm_overlap <ovl> - Sets overlap 10454b9ad928SBarry Smith 10464b9ad928SBarry Smith Notes: 10474b9ad928SBarry Smith By default the ASM preconditioner uses 1 block per processor. To use 10484b9ad928SBarry Smith multiple blocks per perocessor, see PCASMSetTotalSubdomains() and 10494b9ad928SBarry Smith PCASMSetLocalSubdomains() (and the option -pc_asm_blocks <blks>). 10504b9ad928SBarry Smith 10514b9ad928SBarry Smith The overlap defaults to 1, so if one desires that no additional 10524b9ad928SBarry Smith overlap be computed beyond what may have been set with a call to 10534b9ad928SBarry Smith PCASMSetTotalSubdomains() or PCASMSetLocalSubdomains(), then ovl 10544b9ad928SBarry Smith must be set to be 0. In particular, if one does not explicitly set 10554b9ad928SBarry Smith the subdomains an application code, then all overlap would be computed 10564b9ad928SBarry Smith internally by PETSc, and using an overlap of 0 would result in an ASM 10574b9ad928SBarry Smith variant that is equivalent to the block Jacobi preconditioner. 10584b9ad928SBarry Smith 1059f1ee410cSBarry Smith The default algorithm used by PETSc to increase overlap is fast, but not scalable, 1060f1ee410cSBarry Smith use the option -mat_increase_overlap_scalable when the problem and number of processes is large. 1061f1ee410cSBarry Smith 10624b9ad928SBarry Smith Note that one can define initial index sets with any overlap via 1063f1ee410cSBarry Smith PCASMSetLocalSubdomains(); the routine 10644b9ad928SBarry Smith PCASMSetOverlap() merely allows PETSc to extend that overlap further 10654b9ad928SBarry Smith if desired. 10664b9ad928SBarry Smith 10674b9ad928SBarry Smith Level: intermediate 10684b9ad928SBarry Smith 1069db781477SPatrick Sanan .seealso: `PCASMSetTotalSubdomains()`, `PCASMSetLocalSubdomains()`, `PCASMGetSubKSP()`, 1070db781477SPatrick Sanan `PCASMCreateSubdomains2D()`, `PCASMGetLocalSubdomains()`, `MatIncreaseOverlap()` 10714b9ad928SBarry Smith @*/ 10727087cfbeSBarry Smith PetscErrorCode PCASMSetOverlap(PC pc,PetscInt ovl) 10734b9ad928SBarry Smith { 10744b9ad928SBarry Smith PetscFunctionBegin; 10750700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1076c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,ovl,2); 1077cac4c232SBarry Smith PetscTryMethod(pc,"PCASMSetOverlap_C",(PC,PetscInt),(pc,ovl)); 10784b9ad928SBarry Smith PetscFunctionReturn(0); 10794b9ad928SBarry Smith } 10804b9ad928SBarry Smith 10814b9ad928SBarry Smith /*@ 10824b9ad928SBarry Smith PCASMSetType - Sets the type of restriction and interpolation used 10834b9ad928SBarry Smith for local problems in the additive Schwarz method. 10844b9ad928SBarry Smith 1085d083f849SBarry Smith Logically Collective on pc 10864b9ad928SBarry Smith 10874b9ad928SBarry Smith Input Parameters: 10884b9ad928SBarry Smith + pc - the preconditioner context 10894b9ad928SBarry Smith - type - variant of ASM, one of 10904b9ad928SBarry Smith .vb 10914b9ad928SBarry Smith PC_ASM_BASIC - full interpolation and restriction 1092*82b5ce2aSStefano Zampini PC_ASM_RESTRICT - full restriction, local processor interpolation (default) 10934b9ad928SBarry Smith PC_ASM_INTERPOLATE - full interpolation, local processor restriction 10944b9ad928SBarry Smith PC_ASM_NONE - local processor restriction and interpolation 10954b9ad928SBarry Smith .ve 10964b9ad928SBarry Smith 10974b9ad928SBarry Smith Options Database Key: 10984b9ad928SBarry Smith . -pc_asm_type [basic,restrict,interpolate,none] - Sets ASM type 10994b9ad928SBarry Smith 110095452b02SPatrick Sanan Notes: 110195452b02SPatrick Sanan if the is_local arguments are passed to PCASMSetLocalSubdomains() then they are used when PC_ASM_RESTRICT has been selected 1102f1ee410cSBarry Smith to limit the local processor interpolation 1103f1ee410cSBarry Smith 11044b9ad928SBarry Smith Level: intermediate 11054b9ad928SBarry Smith 1106db781477SPatrick Sanan .seealso: `PCASMSetTotalSubdomains()`, `PCASMSetTotalSubdomains()`, `PCASMGetSubKSP()`, 1107db781477SPatrick Sanan `PCASMCreateSubdomains2D()`, `PCASMType`, `PCASMSetLocalType()`, `PCASMGetLocalType()` 11084b9ad928SBarry Smith @*/ 11097087cfbeSBarry Smith PetscErrorCode PCASMSetType(PC pc,PCASMType type) 11104b9ad928SBarry Smith { 11114b9ad928SBarry Smith PetscFunctionBegin; 11120700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1113c5eb9154SBarry Smith PetscValidLogicalCollectiveEnum(pc,type,2); 1114cac4c232SBarry Smith PetscTryMethod(pc,"PCASMSetType_C",(PC,PCASMType),(pc,type)); 11154b9ad928SBarry Smith PetscFunctionReturn(0); 11164b9ad928SBarry Smith } 11174b9ad928SBarry Smith 1118c60c7ad4SBarry Smith /*@ 1119c60c7ad4SBarry Smith PCASMGetType - Gets the type of restriction and interpolation used 1120c60c7ad4SBarry Smith for local problems in the additive Schwarz method. 1121c60c7ad4SBarry Smith 1122d083f849SBarry Smith Logically Collective on pc 1123c60c7ad4SBarry Smith 1124c60c7ad4SBarry Smith Input Parameter: 1125c60c7ad4SBarry Smith . pc - the preconditioner context 1126c60c7ad4SBarry Smith 1127c60c7ad4SBarry Smith Output Parameter: 1128c60c7ad4SBarry Smith . type - variant of ASM, one of 1129c60c7ad4SBarry Smith 1130c60c7ad4SBarry Smith .vb 1131c60c7ad4SBarry Smith PC_ASM_BASIC - full interpolation and restriction 1132c60c7ad4SBarry Smith PC_ASM_RESTRICT - full restriction, local processor interpolation 1133c60c7ad4SBarry Smith PC_ASM_INTERPOLATE - full interpolation, local processor restriction 1134c60c7ad4SBarry Smith PC_ASM_NONE - local processor restriction and interpolation 1135c60c7ad4SBarry Smith .ve 1136c60c7ad4SBarry Smith 1137c60c7ad4SBarry Smith Options Database Key: 1138c60c7ad4SBarry Smith . -pc_asm_type [basic,restrict,interpolate,none] - Sets ASM type 1139c60c7ad4SBarry Smith 1140c60c7ad4SBarry Smith Level: intermediate 1141c60c7ad4SBarry Smith 1142db781477SPatrick Sanan .seealso: `PCASMSetTotalSubdomains()`, `PCASMSetTotalSubdomains()`, `PCASMGetSubKSP()`, 1143db781477SPatrick Sanan `PCASMCreateSubdomains2D()`, `PCASMType`, `PCASMSetType()`, `PCASMSetLocalType()`, `PCASMGetLocalType()` 1144c60c7ad4SBarry Smith @*/ 1145c60c7ad4SBarry Smith PetscErrorCode PCASMGetType(PC pc,PCASMType *type) 1146c60c7ad4SBarry Smith { 1147c60c7ad4SBarry Smith PetscFunctionBegin; 1148c60c7ad4SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1149cac4c232SBarry Smith PetscUseMethod(pc,"PCASMGetType_C",(PC,PCASMType*),(pc,type)); 1150c60c7ad4SBarry Smith PetscFunctionReturn(0); 1151c60c7ad4SBarry Smith } 1152c60c7ad4SBarry Smith 115312cd4985SMatthew G. Knepley /*@ 115412cd4985SMatthew G. Knepley PCASMSetLocalType - Sets the type of composition used for local problems in the additive Schwarz method. 115512cd4985SMatthew G. Knepley 1156d083f849SBarry Smith Logically Collective on pc 115712cd4985SMatthew G. Knepley 115812cd4985SMatthew G. Knepley Input Parameters: 115912cd4985SMatthew G. Knepley + pc - the preconditioner context 116012cd4985SMatthew G. Knepley - type - type of composition, one of 116112cd4985SMatthew G. Knepley .vb 116212cd4985SMatthew G. Knepley PC_COMPOSITE_ADDITIVE - local additive combination 116312cd4985SMatthew G. Knepley PC_COMPOSITE_MULTIPLICATIVE - local multiplicative combination 116412cd4985SMatthew G. Knepley .ve 116512cd4985SMatthew G. Knepley 116612cd4985SMatthew G. Knepley Options Database Key: 116712cd4985SMatthew G. Knepley . -pc_asm_local_type [additive,multiplicative] - Sets local solver composition type 116812cd4985SMatthew G. Knepley 116912cd4985SMatthew G. Knepley Level: intermediate 117012cd4985SMatthew G. Knepley 1171db781477SPatrick Sanan .seealso: `PCASMSetType()`, `PCASMGetType()`, `PCASMGetLocalType()`, `PCASM`, `PCASMType`, `PCASMSetType()`, `PCASMGetType()`, `PCCompositeType` 117212cd4985SMatthew G. Knepley @*/ 117312cd4985SMatthew G. Knepley PetscErrorCode PCASMSetLocalType(PC pc, PCCompositeType type) 117412cd4985SMatthew G. Knepley { 117512cd4985SMatthew G. Knepley PetscFunctionBegin; 117612cd4985SMatthew G. Knepley PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 117712cd4985SMatthew G. Knepley PetscValidLogicalCollectiveEnum(pc, type, 2); 1178cac4c232SBarry Smith PetscTryMethod(pc, "PCASMSetLocalType_C", (PC, PCCompositeType), (pc, type)); 117912cd4985SMatthew G. Knepley PetscFunctionReturn(0); 118012cd4985SMatthew G. Knepley } 118112cd4985SMatthew G. Knepley 118212cd4985SMatthew G. Knepley /*@ 118312cd4985SMatthew G. Knepley PCASMGetLocalType - Gets the type of composition used for local problems in the additive Schwarz method. 118412cd4985SMatthew G. Knepley 1185d083f849SBarry Smith Logically Collective on pc 118612cd4985SMatthew G. Knepley 118712cd4985SMatthew G. Knepley Input Parameter: 118812cd4985SMatthew G. Knepley . pc - the preconditioner context 118912cd4985SMatthew G. Knepley 119012cd4985SMatthew G. Knepley Output Parameter: 119112cd4985SMatthew G. Knepley . type - type of composition, one of 119212cd4985SMatthew G. Knepley .vb 119312cd4985SMatthew G. Knepley PC_COMPOSITE_ADDITIVE - local additive combination 119412cd4985SMatthew G. Knepley PC_COMPOSITE_MULTIPLICATIVE - local multiplicative combination 119512cd4985SMatthew G. Knepley .ve 119612cd4985SMatthew G. Knepley 119712cd4985SMatthew G. Knepley Options Database Key: 119812cd4985SMatthew G. Knepley . -pc_asm_local_type [additive,multiplicative] - Sets local solver composition type 119912cd4985SMatthew G. Knepley 120012cd4985SMatthew G. Knepley Level: intermediate 120112cd4985SMatthew G. Knepley 1202db781477SPatrick Sanan .seealso: `PCASMSetType()`, `PCASMGetType()`, `PCASMSetLocalType()`, `PCASMCreate()`, `PCASMType`, `PCASMSetType()`, `PCASMGetType()`, `PCCompositeType` 120312cd4985SMatthew G. Knepley @*/ 120412cd4985SMatthew G. Knepley PetscErrorCode PCASMGetLocalType(PC pc, PCCompositeType *type) 120512cd4985SMatthew G. Knepley { 120612cd4985SMatthew G. Knepley PetscFunctionBegin; 120712cd4985SMatthew G. Knepley PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 120812cd4985SMatthew G. Knepley PetscValidPointer(type, 2); 1209cac4c232SBarry Smith PetscUseMethod(pc, "PCASMGetLocalType_C", (PC, PCCompositeType *), (pc, type)); 121012cd4985SMatthew G. Knepley PetscFunctionReturn(0); 121112cd4985SMatthew G. Knepley } 121212cd4985SMatthew G. Knepley 12136ed231c7SMatthew Knepley /*@ 12146ed231c7SMatthew Knepley PCASMSetSortIndices - Determines whether subdomain indices are sorted. 12156ed231c7SMatthew Knepley 1216d083f849SBarry Smith Logically Collective on pc 12176ed231c7SMatthew Knepley 12186ed231c7SMatthew Knepley Input Parameters: 12196ed231c7SMatthew Knepley + pc - the preconditioner context 12206ed231c7SMatthew Knepley - doSort - sort the subdomain indices 12216ed231c7SMatthew Knepley 12226ed231c7SMatthew Knepley Level: intermediate 12236ed231c7SMatthew Knepley 1224db781477SPatrick Sanan .seealso: `PCASMSetLocalSubdomains()`, `PCASMSetTotalSubdomains()`, `PCASMGetSubKSP()`, 1225db781477SPatrick Sanan `PCASMCreateSubdomains2D()` 12266ed231c7SMatthew Knepley @*/ 12277087cfbeSBarry Smith PetscErrorCode PCASMSetSortIndices(PC pc,PetscBool doSort) 12286ed231c7SMatthew Knepley { 12296ed231c7SMatthew Knepley PetscFunctionBegin; 12300700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1231acfcf0e5SJed Brown PetscValidLogicalCollectiveBool(pc,doSort,2); 1232cac4c232SBarry Smith PetscTryMethod(pc,"PCASMSetSortIndices_C",(PC,PetscBool),(pc,doSort)); 12336ed231c7SMatthew Knepley PetscFunctionReturn(0); 12346ed231c7SMatthew Knepley } 12356ed231c7SMatthew Knepley 12364b9ad928SBarry Smith /*@C 12374b9ad928SBarry Smith PCASMGetSubKSP - Gets the local KSP contexts for all blocks on 12384b9ad928SBarry Smith this processor. 12394b9ad928SBarry Smith 1240d083f849SBarry Smith Collective on pc iff first_local is requested 12414b9ad928SBarry Smith 12424b9ad928SBarry Smith Input Parameter: 12434b9ad928SBarry Smith . pc - the preconditioner context 12444b9ad928SBarry Smith 12454b9ad928SBarry Smith Output Parameters: 12460298fd71SBarry Smith + n_local - the number of blocks on this processor or NULL 12470298fd71SBarry Smith . first_local - the global number of the first block on this processor or NULL, 12480298fd71SBarry Smith all processors must request or all must pass NULL 12494b9ad928SBarry Smith - ksp - the array of KSP contexts 12504b9ad928SBarry Smith 12514b9ad928SBarry Smith Note: 1252d29017ddSJed Brown After PCASMGetSubKSP() the array of KSPes is not to be freed. 12534b9ad928SBarry Smith 12544b9ad928SBarry Smith You must call KSPSetUp() before calling PCASMGetSubKSP(). 12554b9ad928SBarry Smith 1256d29017ddSJed Brown Fortran note: 12572bf68e3eSBarry 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. 1258d29017ddSJed Brown 12594b9ad928SBarry Smith Level: advanced 12604b9ad928SBarry Smith 1261db781477SPatrick Sanan .seealso: `PCASMSetTotalSubdomains()`, `PCASMSetTotalSubdomains()`, `PCASMSetOverlap()`, 1262db781477SPatrick Sanan `PCASMCreateSubdomains2D()`, 12634b9ad928SBarry Smith @*/ 12647087cfbeSBarry Smith PetscErrorCode PCASMGetSubKSP(PC pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp[]) 12654b9ad928SBarry Smith { 12664b9ad928SBarry Smith PetscFunctionBegin; 12670700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1268cac4c232SBarry Smith PetscUseMethod(pc,"PCASMGetSubKSP_C",(PC,PetscInt*,PetscInt*,KSP **),(pc,n_local,first_local,ksp)); 12694b9ad928SBarry Smith PetscFunctionReturn(0); 12704b9ad928SBarry Smith } 12714b9ad928SBarry Smith 12724b9ad928SBarry Smith /* -------------------------------------------------------------------------------------*/ 12734b9ad928SBarry Smith /*MC 12743b09bd56SBarry Smith PCASM - Use the (restricted) additive Schwarz method, each block is (approximately) solved with 12754b9ad928SBarry Smith its own KSP object. 12764b9ad928SBarry Smith 12774b9ad928SBarry Smith Options Database Keys: 127849517cdeSBarry Smith + -pc_asm_blocks <blks> - Sets total blocks 12794b9ad928SBarry Smith . -pc_asm_overlap <ovl> - Sets overlap 1280f1ee410cSBarry Smith . -pc_asm_type [basic,restrict,interpolate,none] - Sets ASM type, default is restrict 1281f1ee410cSBarry Smith - -pc_asm_local_type [additive, multiplicative] - Sets ASM type, default is additive 12824b9ad928SBarry Smith 12833b09bd56SBarry Smith IMPORTANT: If you run with, for example, 3 blocks on 1 processor or 3 blocks on 3 processors you 12843b09bd56SBarry Smith will get a different convergence rate due to the default option of -pc_asm_type restrict. Use 12853b09bd56SBarry Smith -pc_asm_type basic to use the standard ASM. 12863b09bd56SBarry Smith 128795452b02SPatrick Sanan Notes: 128895452b02SPatrick Sanan Each processor can have one or more blocks, but a block cannot be shared by more 1289f1ee410cSBarry Smith than one processor. Use PCGASM for subdomains shared by multiple processes. Defaults to one block per processor. 12904b9ad928SBarry Smith 12913b09bd56SBarry Smith To set options on the solvers for each block append -sub_ to all the KSP, and PC 1292d7ee0231SBarry Smith options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly 12934b9ad928SBarry Smith 1294a8c7a070SBarry Smith To set the options on the solvers separate for each block call PCASMGetSubKSP() 12954b9ad928SBarry Smith and set the options directly on the resulting KSP object (you can access its PC 12964b9ad928SBarry Smith with KSPGetPC()) 12974b9ad928SBarry Smith 12984b9ad928SBarry Smith Level: beginner 12994b9ad928SBarry Smith 1300c582cd25SBarry Smith References: 1301606c0280SSatish Balay + * - M Dryja, OB Widlund, An additive variant of the Schwarz alternating method for the case of many subregions 130296a0c994SBarry Smith Courant Institute, New York University Technical report 1303606c0280SSatish Balay - * - Barry Smith, Petter Bjorstad, and William Gropp, Domain Decompositions: Parallel Multilevel Methods for Elliptic Partial Differential Equations, 130496a0c994SBarry Smith Cambridge University Press. 1305c582cd25SBarry Smith 1306db781477SPatrick Sanan .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, 1307db781477SPatrick Sanan `PCBJACOBI`, `PCASMGetSubKSP()`, `PCASMSetLocalSubdomains()`, `PCASMType`, `PCASMGetType()`, `PCASMSetLocalType()`, `PCASMGetLocalType()` 1308db781477SPatrick Sanan `PCASMSetTotalSubdomains()`, `PCSetModifySubMatrices()`, `PCASMSetOverlap()`, `PCASMSetType()`, `PCCompositeType` 1309e09e08ccSBarry Smith 13104b9ad928SBarry Smith M*/ 13114b9ad928SBarry Smith 13128cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_ASM(PC pc) 13134b9ad928SBarry Smith { 13144b9ad928SBarry Smith PC_ASM *osm; 13154b9ad928SBarry Smith 13164b9ad928SBarry Smith PetscFunctionBegin; 13179566063dSJacob Faibussowitsch PetscCall(PetscNewLog(pc,&osm)); 13182fa5cd67SKarl Rupp 13194b9ad928SBarry Smith osm->n = PETSC_DECIDE; 13204b9ad928SBarry Smith osm->n_local = 0; 13212b837212SDmitry Karpeev osm->n_local_true = PETSC_DECIDE; 13224b9ad928SBarry Smith osm->overlap = 1; 13230a545947SLisandro Dalcin osm->ksp = NULL; 13240a545947SLisandro Dalcin osm->restriction = NULL; 13250a545947SLisandro Dalcin osm->lprolongation = NULL; 13260a545947SLisandro Dalcin osm->lrestriction = NULL; 13270a545947SLisandro Dalcin osm->x = NULL; 13280a545947SLisandro Dalcin osm->y = NULL; 13290a545947SLisandro Dalcin osm->is = NULL; 13300a545947SLisandro Dalcin osm->is_local = NULL; 13310a545947SLisandro Dalcin osm->mat = NULL; 13320a545947SLisandro Dalcin osm->pmat = NULL; 13334b9ad928SBarry Smith osm->type = PC_ASM_RESTRICT; 133412cd4985SMatthew G. Knepley osm->loctype = PC_COMPOSITE_ADDITIVE; 13356ed231c7SMatthew Knepley osm->sort_indices = PETSC_TRUE; 1336d709ea83SDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 133780ec0b7dSPatrick Sanan osm->sub_mat_type = NULL; 13384b9ad928SBarry Smith 133992bb6962SLisandro Dalcin pc->data = (void*)osm; 13404b9ad928SBarry Smith pc->ops->apply = PCApply_ASM; 134148e38a8aSPierre Jolivet pc->ops->matapply = PCMatApply_ASM; 13424b9ad928SBarry Smith pc->ops->applytranspose = PCApplyTranspose_ASM; 13434b9ad928SBarry Smith pc->ops->setup = PCSetUp_ASM; 1344e91c6855SBarry Smith pc->ops->reset = PCReset_ASM; 13454b9ad928SBarry Smith pc->ops->destroy = PCDestroy_ASM; 13464b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_ASM; 13474b9ad928SBarry Smith pc->ops->setuponblocks = PCSetUpOnBlocks_ASM; 13484b9ad928SBarry Smith pc->ops->view = PCView_ASM; 13490a545947SLisandro Dalcin pc->ops->applyrichardson = NULL; 13504b9ad928SBarry Smith 13519566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCASMSetLocalSubdomains_C",PCASMSetLocalSubdomains_ASM)); 13529566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCASMSetTotalSubdomains_C",PCASMSetTotalSubdomains_ASM)); 13539566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCASMSetOverlap_C",PCASMSetOverlap_ASM)); 13549566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCASMSetType_C",PCASMSetType_ASM)); 13559566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCASMGetType_C",PCASMGetType_ASM)); 13569566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCASMSetLocalType_C",PCASMSetLocalType_ASM)); 13579566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCASMGetLocalType_C",PCASMGetLocalType_ASM)); 13589566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCASMSetSortIndices_C",PCASMSetSortIndices_ASM)); 13599566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCASMGetSubKSP_C",PCASMGetSubKSP_ASM)); 13609566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCASMGetSubMatType_C",PCASMGetSubMatType_ASM)); 13619566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCASMSetSubMatType_C",PCASMSetSubMatType_ASM)); 13624b9ad928SBarry Smith PetscFunctionReturn(0); 13634b9ad928SBarry Smith } 13644b9ad928SBarry Smith 136592bb6962SLisandro Dalcin /*@C 136692bb6962SLisandro Dalcin PCASMCreateSubdomains - Creates the index sets for the overlapping Schwarz 13673b8d980fSPierre Jolivet preconditioner for any problem on a general grid. 136892bb6962SLisandro Dalcin 136992bb6962SLisandro Dalcin Collective 137092bb6962SLisandro Dalcin 137192bb6962SLisandro Dalcin Input Parameters: 137292bb6962SLisandro Dalcin + A - The global matrix operator 137392bb6962SLisandro Dalcin - n - the number of local blocks 137492bb6962SLisandro Dalcin 137592bb6962SLisandro Dalcin Output Parameters: 137692bb6962SLisandro Dalcin . outis - the array of index sets defining the subdomains 137792bb6962SLisandro Dalcin 137892bb6962SLisandro Dalcin Level: advanced 137992bb6962SLisandro Dalcin 13807d6bfa3bSBarry Smith Note: this generates nonoverlapping subdomains; the PCASM will generate the overlap 13817d6bfa3bSBarry Smith from these if you use PCASMSetLocalSubdomains() 13827d6bfa3bSBarry Smith 13837d6bfa3bSBarry Smith In the Fortran version you must provide the array outis[] already allocated of length n. 13847d6bfa3bSBarry Smith 1385db781477SPatrick Sanan .seealso: `PCASMSetLocalSubdomains()`, `PCASMDestroySubdomains()` 138692bb6962SLisandro Dalcin @*/ 13877087cfbeSBarry Smith PetscErrorCode PCASMCreateSubdomains(Mat A, PetscInt n, IS* outis[]) 138892bb6962SLisandro Dalcin { 138992bb6962SLisandro Dalcin MatPartitioning mpart; 139092bb6962SLisandro Dalcin const char *prefix; 139192bb6962SLisandro Dalcin PetscInt i,j,rstart,rend,bs; 1392976e8c5aSLisandro Dalcin PetscBool hasop, isbaij = PETSC_FALSE,foundpart = PETSC_FALSE; 13930298fd71SBarry Smith Mat Ad = NULL, adj; 139492bb6962SLisandro Dalcin IS ispart,isnumb,*is; 139592bb6962SLisandro Dalcin 139692bb6962SLisandro Dalcin PetscFunctionBegin; 13970700a824SBarry Smith PetscValidHeaderSpecific(A,MAT_CLASSID,1); 139892bb6962SLisandro Dalcin PetscValidPointer(outis,3); 139963a3b9bcSJacob Faibussowitsch PetscCheck(n >= 1,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of local blocks must be > 0, n = %" PetscInt_FMT,n); 140092bb6962SLisandro Dalcin 140192bb6962SLisandro Dalcin /* Get prefix, row distribution, and block size */ 14029566063dSJacob Faibussowitsch PetscCall(MatGetOptionsPrefix(A,&prefix)); 14039566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A,&rstart,&rend)); 14049566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(A,&bs)); 140563a3b9bcSJacob 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); 140665e19b50SBarry Smith 140792bb6962SLisandro Dalcin /* Get diagonal block from matrix if possible */ 14089566063dSJacob Faibussowitsch PetscCall(MatHasOperation(A,MATOP_GET_DIAGONAL_BLOCK,&hasop)); 1409976e8c5aSLisandro Dalcin if (hasop) { 14109566063dSJacob Faibussowitsch PetscCall(MatGetDiagonalBlock(A,&Ad)); 141192bb6962SLisandro Dalcin } 141292bb6962SLisandro Dalcin if (Ad) { 14139566063dSJacob Faibussowitsch PetscCall(PetscObjectBaseTypeCompare((PetscObject)Ad,MATSEQBAIJ,&isbaij)); 14149566063dSJacob Faibussowitsch if (!isbaij) PetscCall(PetscObjectBaseTypeCompare((PetscObject)Ad,MATSEQSBAIJ,&isbaij)); 141592bb6962SLisandro Dalcin } 141692bb6962SLisandro Dalcin if (Ad && n > 1) { 1417ace3abfcSBarry Smith PetscBool match,done; 141892bb6962SLisandro Dalcin /* Try to setup a good matrix partitioning if available */ 14199566063dSJacob Faibussowitsch PetscCall(MatPartitioningCreate(PETSC_COMM_SELF,&mpart)); 14209566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix)); 14219566063dSJacob Faibussowitsch PetscCall(MatPartitioningSetFromOptions(mpart)); 14229566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGCURRENT,&match)); 142392bb6962SLisandro Dalcin if (!match) { 14249566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGSQUARE,&match)); 142592bb6962SLisandro Dalcin } 142692bb6962SLisandro Dalcin if (!match) { /* assume a "good" partitioner is available */ 14271a83f524SJed Brown PetscInt na; 14281a83f524SJed Brown const PetscInt *ia,*ja; 14299566063dSJacob Faibussowitsch PetscCall(MatGetRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done)); 143092bb6962SLisandro Dalcin if (done) { 143192bb6962SLisandro Dalcin /* Build adjacency matrix by hand. Unfortunately a call to 143292bb6962SLisandro Dalcin MatConvert(Ad,MATMPIADJ,MAT_INITIAL_MATRIX,&adj) will 143392bb6962SLisandro Dalcin remove the block-aij structure and we cannot expect 143492bb6962SLisandro Dalcin MatPartitioning to split vertices as we need */ 14350a545947SLisandro Dalcin PetscInt i,j,len,nnz,cnt,*iia=NULL,*jja=NULL; 14361a83f524SJed Brown const PetscInt *row; 143792bb6962SLisandro Dalcin nnz = 0; 143892bb6962SLisandro Dalcin for (i=0; i<na; i++) { /* count number of nonzeros */ 143992bb6962SLisandro Dalcin len = ia[i+1] - ia[i]; 144092bb6962SLisandro Dalcin row = ja + ia[i]; 144192bb6962SLisandro Dalcin for (j=0; j<len; j++) { 144292bb6962SLisandro Dalcin if (row[j] == i) { /* don't count diagonal */ 144392bb6962SLisandro Dalcin len--; break; 144492bb6962SLisandro Dalcin } 144592bb6962SLisandro Dalcin } 144692bb6962SLisandro Dalcin nnz += len; 144792bb6962SLisandro Dalcin } 14489566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(na+1,&iia)); 14499566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nnz,&jja)); 145092bb6962SLisandro Dalcin nnz = 0; 145192bb6962SLisandro Dalcin iia[0] = 0; 145292bb6962SLisandro Dalcin for (i=0; i<na; i++) { /* fill adjacency */ 145392bb6962SLisandro Dalcin cnt = 0; 145492bb6962SLisandro Dalcin len = ia[i+1] - ia[i]; 145592bb6962SLisandro Dalcin row = ja + ia[i]; 145692bb6962SLisandro Dalcin for (j=0; j<len; j++) { 145792bb6962SLisandro Dalcin if (row[j] != i) { /* if not diagonal */ 145892bb6962SLisandro Dalcin jja[nnz+cnt++] = row[j]; 145992bb6962SLisandro Dalcin } 146092bb6962SLisandro Dalcin } 146192bb6962SLisandro Dalcin nnz += cnt; 146292bb6962SLisandro Dalcin iia[i+1] = nnz; 146392bb6962SLisandro Dalcin } 146492bb6962SLisandro Dalcin /* Partitioning of the adjacency matrix */ 14659566063dSJacob Faibussowitsch PetscCall(MatCreateMPIAdj(PETSC_COMM_SELF,na,na,iia,jja,NULL,&adj)); 14669566063dSJacob Faibussowitsch PetscCall(MatPartitioningSetAdjacency(mpart,adj)); 14679566063dSJacob Faibussowitsch PetscCall(MatPartitioningSetNParts(mpart,n)); 14689566063dSJacob Faibussowitsch PetscCall(MatPartitioningApply(mpart,&ispart)); 14699566063dSJacob Faibussowitsch PetscCall(ISPartitioningToNumbering(ispart,&isnumb)); 14709566063dSJacob Faibussowitsch PetscCall(MatDestroy(&adj)); 147192bb6962SLisandro Dalcin foundpart = PETSC_TRUE; 147292bb6962SLisandro Dalcin } 14739566063dSJacob Faibussowitsch PetscCall(MatRestoreRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done)); 147492bb6962SLisandro Dalcin } 14759566063dSJacob Faibussowitsch PetscCall(MatPartitioningDestroy(&mpart)); 147692bb6962SLisandro Dalcin } 147792bb6962SLisandro Dalcin 14789566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n,&is)); 147992bb6962SLisandro Dalcin *outis = is; 148092bb6962SLisandro Dalcin 148192bb6962SLisandro Dalcin if (!foundpart) { 148292bb6962SLisandro Dalcin 148392bb6962SLisandro Dalcin /* Partitioning by contiguous chunks of rows */ 148492bb6962SLisandro Dalcin 148592bb6962SLisandro Dalcin PetscInt mbs = (rend-rstart)/bs; 148692bb6962SLisandro Dalcin PetscInt start = rstart; 148792bb6962SLisandro Dalcin for (i=0; i<n; i++) { 148892bb6962SLisandro Dalcin PetscInt count = (mbs/n + ((mbs % n) > i)) * bs; 14899566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF,count,start,1,&is[i])); 149092bb6962SLisandro Dalcin start += count; 149192bb6962SLisandro Dalcin } 149292bb6962SLisandro Dalcin 149392bb6962SLisandro Dalcin } else { 149492bb6962SLisandro Dalcin 149592bb6962SLisandro Dalcin /* Partitioning by adjacency of diagonal block */ 149692bb6962SLisandro Dalcin 149792bb6962SLisandro Dalcin const PetscInt *numbering; 149892bb6962SLisandro Dalcin PetscInt *count,nidx,*indices,*newidx,start=0; 149992bb6962SLisandro Dalcin /* Get node count in each partition */ 15009566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n,&count)); 15019566063dSJacob Faibussowitsch PetscCall(ISPartitioningCount(ispart,n,count)); 150292bb6962SLisandro Dalcin if (isbaij && bs > 1) { /* adjust for the block-aij case */ 150392bb6962SLisandro Dalcin for (i=0; i<n; i++) count[i] *= bs; 150492bb6962SLisandro Dalcin } 150592bb6962SLisandro Dalcin /* Build indices from node numbering */ 15069566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(isnumb,&nidx)); 15079566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nidx,&indices)); 150892bb6962SLisandro Dalcin for (i=0; i<nidx; i++) indices[i] = i; /* needs to be initialized */ 15099566063dSJacob Faibussowitsch PetscCall(ISGetIndices(isnumb,&numbering)); 15109566063dSJacob Faibussowitsch PetscCall(PetscSortIntWithPermutation(nidx,numbering,indices)); 15119566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(isnumb,&numbering)); 151292bb6962SLisandro Dalcin if (isbaij && bs > 1) { /* adjust for the block-aij case */ 15139566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nidx*bs,&newidx)); 15142fa5cd67SKarl Rupp for (i=0; i<nidx; i++) { 15152fa5cd67SKarl Rupp for (j=0; j<bs; j++) newidx[i*bs+j] = indices[i]*bs + j; 15162fa5cd67SKarl Rupp } 15179566063dSJacob Faibussowitsch PetscCall(PetscFree(indices)); 151892bb6962SLisandro Dalcin nidx *= bs; 151992bb6962SLisandro Dalcin indices = newidx; 152092bb6962SLisandro Dalcin } 152192bb6962SLisandro Dalcin /* Shift to get global indices */ 152292bb6962SLisandro Dalcin for (i=0; i<nidx; i++) indices[i] += rstart; 152392bb6962SLisandro Dalcin 152492bb6962SLisandro Dalcin /* Build the index sets for each block */ 152592bb6962SLisandro Dalcin for (i=0; i<n; i++) { 15269566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF,count[i],&indices[start],PETSC_COPY_VALUES,&is[i])); 15279566063dSJacob Faibussowitsch PetscCall(ISSort(is[i])); 152892bb6962SLisandro Dalcin start += count[i]; 152992bb6962SLisandro Dalcin } 153092bb6962SLisandro Dalcin 15319566063dSJacob Faibussowitsch PetscCall(PetscFree(count)); 15329566063dSJacob Faibussowitsch PetscCall(PetscFree(indices)); 15339566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isnumb)); 15349566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ispart)); 153592bb6962SLisandro Dalcin 153692bb6962SLisandro Dalcin } 153792bb6962SLisandro Dalcin PetscFunctionReturn(0); 153892bb6962SLisandro Dalcin } 153992bb6962SLisandro Dalcin 154092bb6962SLisandro Dalcin /*@C 154192bb6962SLisandro Dalcin PCASMDestroySubdomains - Destroys the index sets created with 154292bb6962SLisandro Dalcin PCASMCreateSubdomains(). Should be called after setting subdomains 154392bb6962SLisandro Dalcin with PCASMSetLocalSubdomains(). 154492bb6962SLisandro Dalcin 154592bb6962SLisandro Dalcin Collective 154692bb6962SLisandro Dalcin 154792bb6962SLisandro Dalcin Input Parameters: 154892bb6962SLisandro Dalcin + n - the number of index sets 15492b691e39SMatthew Knepley . is - the array of index sets 15500298fd71SBarry Smith - is_local - the array of local index sets, can be NULL 155192bb6962SLisandro Dalcin 155292bb6962SLisandro Dalcin Level: advanced 155392bb6962SLisandro Dalcin 1554db781477SPatrick Sanan .seealso: `PCASMCreateSubdomains()`, `PCASMSetLocalSubdomains()` 155592bb6962SLisandro Dalcin @*/ 15567087cfbeSBarry Smith PetscErrorCode PCASMDestroySubdomains(PetscInt n, IS is[], IS is_local[]) 155792bb6962SLisandro Dalcin { 155892bb6962SLisandro Dalcin PetscInt i; 15595fd66863SKarl Rupp 156092bb6962SLisandro Dalcin PetscFunctionBegin; 1561a35b7b57SDmitry Karpeev if (n <= 0) PetscFunctionReturn(0); 1562a35b7b57SDmitry Karpeev if (is) { 156392bb6962SLisandro Dalcin PetscValidPointer(is,2); 15649566063dSJacob Faibussowitsch for (i=0; i<n; i++) PetscCall(ISDestroy(&is[i])); 15659566063dSJacob Faibussowitsch PetscCall(PetscFree(is)); 1566a35b7b57SDmitry Karpeev } 15672b691e39SMatthew Knepley if (is_local) { 15682b691e39SMatthew Knepley PetscValidPointer(is_local,3); 15699566063dSJacob Faibussowitsch for (i=0; i<n; i++) PetscCall(ISDestroy(&is_local[i])); 15709566063dSJacob Faibussowitsch PetscCall(PetscFree(is_local)); 15712b691e39SMatthew Knepley } 157292bb6962SLisandro Dalcin PetscFunctionReturn(0); 157392bb6962SLisandro Dalcin } 157492bb6962SLisandro Dalcin 15754b9ad928SBarry Smith /*@ 15764b9ad928SBarry Smith PCASMCreateSubdomains2D - Creates the index sets for the overlapping Schwarz 15774b9ad928SBarry Smith preconditioner for a two-dimensional problem on a regular grid. 15784b9ad928SBarry Smith 15794b9ad928SBarry Smith Not Collective 15804b9ad928SBarry Smith 15814b9ad928SBarry Smith Input Parameters: 15826b867d5aSJose E. Roman + m - the number of mesh points in the x direction 15836b867d5aSJose E. Roman . n - the number of mesh points in the y direction 15846b867d5aSJose E. Roman . M - the number of subdomains in the x direction 15856b867d5aSJose E. Roman . N - the number of subdomains in the y direction 15864b9ad928SBarry Smith . dof - degrees of freedom per node 15874b9ad928SBarry Smith - overlap - overlap in mesh lines 15884b9ad928SBarry Smith 15894b9ad928SBarry Smith Output Parameters: 15904b9ad928SBarry Smith + Nsub - the number of subdomains created 15913d03bb48SJed Brown . is - array of index sets defining overlapping (if overlap > 0) subdomains 15923d03bb48SJed Brown - is_local - array of index sets defining non-overlapping subdomains 15934b9ad928SBarry Smith 15944b9ad928SBarry Smith Note: 15954b9ad928SBarry Smith Presently PCAMSCreateSubdomains2d() is valid only for sequential 15964b9ad928SBarry Smith preconditioners. More general related routines are 15974b9ad928SBarry Smith PCASMSetTotalSubdomains() and PCASMSetLocalSubdomains(). 15984b9ad928SBarry Smith 15994b9ad928SBarry Smith Level: advanced 16004b9ad928SBarry Smith 1601db781477SPatrick Sanan .seealso: `PCASMSetTotalSubdomains()`, `PCASMSetLocalSubdomains()`, `PCASMGetSubKSP()`, 1602db781477SPatrick Sanan `PCASMSetOverlap()` 16034b9ad928SBarry Smith @*/ 16047087cfbeSBarry Smith PetscErrorCode PCASMCreateSubdomains2D(PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt dof,PetscInt overlap,PetscInt *Nsub,IS **is,IS **is_local) 16054b9ad928SBarry Smith { 16063d03bb48SJed Brown PetscInt i,j,height,width,ystart,xstart,yleft,yright,xleft,xright,loc_outer; 160713f74950SBarry Smith PetscInt nidx,*idx,loc,ii,jj,count; 16084b9ad928SBarry Smith 16094b9ad928SBarry Smith PetscFunctionBegin; 16107827d75bSBarry Smith PetscCheck(dof == 1,PETSC_COMM_SELF,PETSC_ERR_SUP,"dof must be 1"); 16114b9ad928SBarry Smith 16124b9ad928SBarry Smith *Nsub = N*M; 16139566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(*Nsub,is)); 16149566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(*Nsub,is_local)); 16154b9ad928SBarry Smith ystart = 0; 16163d03bb48SJed Brown loc_outer = 0; 16174b9ad928SBarry Smith for (i=0; i<N; i++) { 16184b9ad928SBarry Smith height = n/N + ((n % N) > i); /* height of subdomain */ 16197827d75bSBarry Smith PetscCheck(height >= 2,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many N subdomains for mesh dimension n"); 16204b9ad928SBarry Smith yleft = ystart - overlap; if (yleft < 0) yleft = 0; 16214b9ad928SBarry Smith yright = ystart + height + overlap; if (yright > n) yright = n; 16224b9ad928SBarry Smith xstart = 0; 16234b9ad928SBarry Smith for (j=0; j<M; j++) { 16244b9ad928SBarry Smith width = m/M + ((m % M) > j); /* width of subdomain */ 16257827d75bSBarry Smith PetscCheck(width >= 2,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many M subdomains for mesh dimension m"); 16264b9ad928SBarry Smith xleft = xstart - overlap; if (xleft < 0) xleft = 0; 16274b9ad928SBarry Smith xright = xstart + width + overlap; if (xright > m) xright = m; 16284b9ad928SBarry Smith nidx = (xright - xleft)*(yright - yleft); 16299566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nidx,&idx)); 16304b9ad928SBarry Smith loc = 0; 16314b9ad928SBarry Smith for (ii=yleft; ii<yright; ii++) { 16324b9ad928SBarry Smith count = m*ii + xleft; 16332fa5cd67SKarl Rupp for (jj=xleft; jj<xright; jj++) idx[loc++] = count++; 16344b9ad928SBarry Smith } 16359566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF,nidx,idx,PETSC_COPY_VALUES,(*is)+loc_outer)); 16363d03bb48SJed Brown if (overlap == 0) { 16379566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)(*is)[loc_outer])); 16382fa5cd67SKarl Rupp 16393d03bb48SJed Brown (*is_local)[loc_outer] = (*is)[loc_outer]; 16403d03bb48SJed Brown } else { 16413d03bb48SJed Brown for (loc=0,ii=ystart; ii<ystart+height; ii++) { 16423d03bb48SJed Brown for (jj=xstart; jj<xstart+width; jj++) { 16433d03bb48SJed Brown idx[loc++] = m*ii + jj; 16443d03bb48SJed Brown } 16453d03bb48SJed Brown } 16469566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF,loc,idx,PETSC_COPY_VALUES,*is_local+loc_outer)); 16473d03bb48SJed Brown } 16489566063dSJacob Faibussowitsch PetscCall(PetscFree(idx)); 16494b9ad928SBarry Smith xstart += width; 16503d03bb48SJed Brown loc_outer++; 16514b9ad928SBarry Smith } 16524b9ad928SBarry Smith ystart += height; 16534b9ad928SBarry Smith } 16549566063dSJacob Faibussowitsch for (i=0; i<*Nsub; i++) PetscCall(ISSort((*is)[i])); 16554b9ad928SBarry Smith PetscFunctionReturn(0); 16564b9ad928SBarry Smith } 16574b9ad928SBarry Smith 16584b9ad928SBarry Smith /*@C 16594b9ad928SBarry Smith PCASMGetLocalSubdomains - Gets the local subdomains (for this processor 16604b9ad928SBarry Smith only) for the additive Schwarz preconditioner. 16614b9ad928SBarry Smith 1662ad4df100SBarry Smith Not Collective 16634b9ad928SBarry Smith 16644b9ad928SBarry Smith Input Parameter: 16654b9ad928SBarry Smith . pc - the preconditioner context 16664b9ad928SBarry Smith 16674b9ad928SBarry Smith Output Parameters: 1668f17a6784SPierre Jolivet + n - if requested, the number of subdomains for this processor (default value = 1) 1669f17a6784SPierre Jolivet . is - if requested, the index sets that define the subdomains for this processor 1670f17a6784SPierre Jolivet - is_local - if requested, the index sets that define the local part of the subdomains for this processor (can be NULL) 16714b9ad928SBarry Smith 16724b9ad928SBarry Smith Notes: 16734b9ad928SBarry Smith The IS numbering is in the parallel, global numbering of the vector. 16744b9ad928SBarry Smith 16754b9ad928SBarry Smith Level: advanced 16764b9ad928SBarry Smith 1677db781477SPatrick Sanan .seealso: `PCASMSetTotalSubdomains()`, `PCASMSetOverlap()`, `PCASMGetSubKSP()`, 1678db781477SPatrick Sanan `PCASMCreateSubdomains2D()`, `PCASMSetLocalSubdomains()`, `PCASMGetLocalSubmatrices()` 16794b9ad928SBarry Smith @*/ 16807087cfbeSBarry Smith PetscErrorCode PCASMGetLocalSubdomains(PC pc,PetscInt *n,IS *is[],IS *is_local[]) 16814b9ad928SBarry Smith { 16822a808120SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 1683ace3abfcSBarry Smith PetscBool match; 16844b9ad928SBarry Smith 16854b9ad928SBarry Smith PetscFunctionBegin; 16860700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1687f17a6784SPierre Jolivet if (n) PetscValidIntPointer(n,2); 168892bb6962SLisandro Dalcin if (is) PetscValidPointer(is,3); 1689f17a6784SPierre Jolivet if (is_local) PetscValidPointer(is_local,4); 16909566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc,PCASM,&match)); 169128b400f6SJacob Faibussowitsch PetscCheck(match,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"PC is not a PCASM"); 16924b9ad928SBarry Smith if (n) *n = osm->n_local_true; 16934b9ad928SBarry Smith if (is) *is = osm->is; 16942b691e39SMatthew Knepley if (is_local) *is_local = osm->is_local; 16954b9ad928SBarry Smith PetscFunctionReturn(0); 16964b9ad928SBarry Smith } 16974b9ad928SBarry Smith 16984b9ad928SBarry Smith /*@C 16994b9ad928SBarry Smith PCASMGetLocalSubmatrices - Gets the local submatrices (for this processor 17004b9ad928SBarry Smith only) for the additive Schwarz preconditioner. 17014b9ad928SBarry Smith 1702ad4df100SBarry Smith Not Collective 17034b9ad928SBarry Smith 17044b9ad928SBarry Smith Input Parameter: 17054b9ad928SBarry Smith . pc - the preconditioner context 17064b9ad928SBarry Smith 17074b9ad928SBarry Smith Output Parameters: 1708f17a6784SPierre Jolivet + n - if requested, the number of matrices for this processor (default value = 1) 1709f17a6784SPierre Jolivet - mat - if requested, the matrices 17104b9ad928SBarry Smith 17114b9ad928SBarry Smith Level: advanced 17124b9ad928SBarry Smith 171395452b02SPatrick Sanan Notes: 171434a84908Sprj- Call after PCSetUp() (or KSPSetUp()) but before PCApply() and before PCSetUpOnBlocks()) 1715cf739d55SBarry Smith 171634a84908Sprj- Usually one would use PCSetModifySubMatrices() to change the submatrices in building the preconditioner. 1717cf739d55SBarry Smith 1718db781477SPatrick Sanan .seealso: `PCASMSetTotalSubdomains()`, `PCASMSetOverlap()`, `PCASMGetSubKSP()`, 1719db781477SPatrick Sanan `PCASMCreateSubdomains2D()`, `PCASMSetLocalSubdomains()`, `PCASMGetLocalSubdomains()`, `PCSetModifySubMatrices()` 17204b9ad928SBarry Smith @*/ 17217087cfbeSBarry Smith PetscErrorCode PCASMGetLocalSubmatrices(PC pc,PetscInt *n,Mat *mat[]) 17224b9ad928SBarry Smith { 17234b9ad928SBarry Smith PC_ASM *osm; 1724ace3abfcSBarry Smith PetscBool match; 17254b9ad928SBarry Smith 17264b9ad928SBarry Smith PetscFunctionBegin; 17270700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1728f17a6784SPierre Jolivet if (n) PetscValidIntPointer(n,2); 172992bb6962SLisandro Dalcin if (mat) PetscValidPointer(mat,3); 173028b400f6SJacob Faibussowitsch PetscCheck(pc->setupcalled,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call after KSPSetUp() or PCSetUp()."); 17319566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc,PCASM,&match)); 173292bb6962SLisandro Dalcin if (!match) { 173392bb6962SLisandro Dalcin if (n) *n = 0; 17340298fd71SBarry Smith if (mat) *mat = NULL; 173592bb6962SLisandro Dalcin } else { 17364b9ad928SBarry Smith osm = (PC_ASM*)pc->data; 17374b9ad928SBarry Smith if (n) *n = osm->n_local_true; 17384b9ad928SBarry Smith if (mat) *mat = osm->pmat; 173992bb6962SLisandro Dalcin } 17404b9ad928SBarry Smith PetscFunctionReturn(0); 17414b9ad928SBarry Smith } 1742d709ea83SDmitry Karpeev 1743d709ea83SDmitry Karpeev /*@ 1744d709ea83SDmitry Karpeev PCASMSetDMSubdomains - Indicates whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible. 1745f1ee410cSBarry Smith 1746d709ea83SDmitry Karpeev Logically Collective 1747d709ea83SDmitry Karpeev 1748d8d19677SJose E. Roman Input Parameters: 1749d709ea83SDmitry Karpeev + pc - the preconditioner 1750d709ea83SDmitry Karpeev - flg - boolean indicating whether to use subdomains defined by the DM 1751d709ea83SDmitry Karpeev 1752d709ea83SDmitry Karpeev Options Database Key: 1753147403d9SBarry Smith . -pc_asm_dm_subdomains <bool> - use subdomains defined by the DM 1754d709ea83SDmitry Karpeev 1755d709ea83SDmitry Karpeev Level: intermediate 1756d709ea83SDmitry Karpeev 1757d709ea83SDmitry Karpeev Notes: 1758d709ea83SDmitry Karpeev PCASMSetTotalSubdomains() and PCASMSetOverlap() take precedence over PCASMSetDMSubdomains(), 1759d709ea83SDmitry Karpeev so setting either of the first two effectively turns the latter off. 1760d709ea83SDmitry Karpeev 1761db781477SPatrick Sanan .seealso: `PCASMGetDMSubdomains()`, `PCASMSetTotalSubdomains()`, `PCASMSetOverlap()` 1762db781477SPatrick Sanan `PCASMCreateSubdomains2D()`, `PCASMSetLocalSubdomains()`, `PCASMGetLocalSubdomains()` 1763d709ea83SDmitry Karpeev @*/ 1764d709ea83SDmitry Karpeev PetscErrorCode PCASMSetDMSubdomains(PC pc,PetscBool flg) 1765d709ea83SDmitry Karpeev { 1766d709ea83SDmitry Karpeev PC_ASM *osm = (PC_ASM*)pc->data; 1767d709ea83SDmitry Karpeev PetscBool match; 1768d709ea83SDmitry Karpeev 1769d709ea83SDmitry Karpeev PetscFunctionBegin; 1770d709ea83SDmitry Karpeev PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1771d709ea83SDmitry Karpeev PetscValidLogicalCollectiveBool(pc,flg,2); 177228b400f6SJacob Faibussowitsch PetscCheck(!pc->setupcalled,((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Not for a setup PC."); 17739566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc,PCASM,&match)); 1774d709ea83SDmitry Karpeev if (match) { 1775d709ea83SDmitry Karpeev osm->dm_subdomains = flg; 1776d709ea83SDmitry Karpeev } 1777d709ea83SDmitry Karpeev PetscFunctionReturn(0); 1778d709ea83SDmitry Karpeev } 1779d709ea83SDmitry Karpeev 1780d709ea83SDmitry Karpeev /*@ 1781d709ea83SDmitry Karpeev PCASMGetDMSubdomains - Returns flag indicating whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible. 1782d709ea83SDmitry Karpeev Not Collective 1783d709ea83SDmitry Karpeev 1784d709ea83SDmitry Karpeev Input Parameter: 1785d709ea83SDmitry Karpeev . pc - the preconditioner 1786d709ea83SDmitry Karpeev 1787d709ea83SDmitry Karpeev Output Parameter: 1788d709ea83SDmitry Karpeev . flg - boolean indicating whether to use subdomains defined by the DM 1789d709ea83SDmitry Karpeev 1790d709ea83SDmitry Karpeev Level: intermediate 1791d709ea83SDmitry Karpeev 1792db781477SPatrick Sanan .seealso: `PCASMSetDMSubdomains()`, `PCASMSetTotalSubdomains()`, `PCASMSetOverlap()` 1793db781477SPatrick Sanan `PCASMCreateSubdomains2D()`, `PCASMSetLocalSubdomains()`, `PCASMGetLocalSubdomains()` 1794d709ea83SDmitry Karpeev @*/ 1795d709ea83SDmitry Karpeev PetscErrorCode PCASMGetDMSubdomains(PC pc,PetscBool* flg) 1796d709ea83SDmitry Karpeev { 1797d709ea83SDmitry Karpeev PC_ASM *osm = (PC_ASM*)pc->data; 1798d709ea83SDmitry Karpeev PetscBool match; 1799d709ea83SDmitry Karpeev 1800d709ea83SDmitry Karpeev PetscFunctionBegin; 1801d709ea83SDmitry Karpeev PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1802534a8f05SLisandro Dalcin PetscValidBoolPointer(flg,2); 18039566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc,PCASM,&match)); 180456ea09ceSStefano Zampini if (match) *flg = osm->dm_subdomains; 180556ea09ceSStefano Zampini else *flg = PETSC_FALSE; 1806d709ea83SDmitry Karpeev PetscFunctionReturn(0); 1807d709ea83SDmitry Karpeev } 180880ec0b7dSPatrick Sanan 180980ec0b7dSPatrick Sanan /*@ 181080ec0b7dSPatrick Sanan PCASMGetSubMatType - Gets the matrix type used for ASM subsolves, as a string. 181180ec0b7dSPatrick Sanan 181280ec0b7dSPatrick Sanan Not Collective 181380ec0b7dSPatrick Sanan 181480ec0b7dSPatrick Sanan Input Parameter: 181580ec0b7dSPatrick Sanan . pc - the PC 181680ec0b7dSPatrick Sanan 181780ec0b7dSPatrick Sanan Output Parameter: 1818f1ee410cSBarry Smith . -pc_asm_sub_mat_type - name of matrix type 181980ec0b7dSPatrick Sanan 182080ec0b7dSPatrick Sanan Level: advanced 182180ec0b7dSPatrick Sanan 1822db781477SPatrick Sanan .seealso: `PCASMSetSubMatType()`, `PCASM`, `PCSetType()`, `VecSetType()`, `MatType`, `Mat` 182380ec0b7dSPatrick Sanan @*/ 182456ea09ceSStefano Zampini PetscErrorCode PCASMGetSubMatType(PC pc,MatType *sub_mat_type) 182556ea09ceSStefano Zampini { 182656ea09ceSStefano Zampini PetscFunctionBegin; 182756ea09ceSStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1828cac4c232SBarry Smith PetscTryMethod(pc,"PCASMGetSubMatType_C",(PC,MatType*),(pc,sub_mat_type)); 182980ec0b7dSPatrick Sanan PetscFunctionReturn(0); 183080ec0b7dSPatrick Sanan } 183180ec0b7dSPatrick Sanan 183280ec0b7dSPatrick Sanan /*@ 183380ec0b7dSPatrick Sanan PCASMSetSubMatType - Set the type of matrix used for ASM subsolves 183480ec0b7dSPatrick Sanan 183580ec0b7dSPatrick Sanan Collective on Mat 183680ec0b7dSPatrick Sanan 183780ec0b7dSPatrick Sanan Input Parameters: 183880ec0b7dSPatrick Sanan + pc - the PC object 183980ec0b7dSPatrick Sanan - sub_mat_type - matrix type 184080ec0b7dSPatrick Sanan 184180ec0b7dSPatrick Sanan Options Database Key: 184280ec0b7dSPatrick 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. 184380ec0b7dSPatrick Sanan 184480ec0b7dSPatrick Sanan Notes: 184580ec0b7dSPatrick Sanan See "${PETSC_DIR}/include/petscmat.h" for available types 184680ec0b7dSPatrick Sanan 184780ec0b7dSPatrick Sanan Level: advanced 184880ec0b7dSPatrick Sanan 1849db781477SPatrick Sanan .seealso: `PCASMGetSubMatType()`, `PCASM`, `PCSetType()`, `VecSetType()`, `MatType`, `Mat` 185080ec0b7dSPatrick Sanan @*/ 185180ec0b7dSPatrick Sanan PetscErrorCode PCASMSetSubMatType(PC pc,MatType sub_mat_type) 185280ec0b7dSPatrick Sanan { 185356ea09ceSStefano Zampini PetscFunctionBegin; 185456ea09ceSStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1855cac4c232SBarry Smith PetscTryMethod(pc,"PCASMSetSubMatType_C",(PC,MatType),(pc,sub_mat_type)); 185680ec0b7dSPatrick Sanan PetscFunctionReturn(0); 185780ec0b7dSPatrick Sanan } 1858