1dba47a55SKris Buschelman 24b9ad928SBarry Smith /* 34b9ad928SBarry Smith This file defines an additive Schwarz preconditioner for any Mat implementation. 44b9ad928SBarry Smith 54b9ad928SBarry Smith Note that each processor may have any number of subdomains. But in order to 64b9ad928SBarry Smith deal easily with the VecScatter(), we treat each processor as if it has the 74b9ad928SBarry Smith same number of subdomains. 84b9ad928SBarry Smith 94b9ad928SBarry Smith n - total number of true subdomains on all processors 104b9ad928SBarry Smith n_local_true - actual number of subdomains on this processor 114b9ad928SBarry Smith n_local = maximum over all processors of n_local_true 124b9ad928SBarry Smith */ 13af0996ceSBarry Smith #include <petsc/private/pcimpl.h> /*I "petscpc.h" I*/ 141e25c274SJed Brown #include <petscdm.h> 154b9ad928SBarry Smith 164b9ad928SBarry Smith typedef struct { 1713f74950SBarry Smith PetscInt n, n_local, n_local_true; 1813f74950SBarry Smith PetscInt overlap; /* overlap requested by user */ 194b9ad928SBarry Smith KSP *ksp; /* linear solvers for each block */ 203c4dc4c4SMatthew Knepley VecScatter *restriction; /* mapping from global to subregion */ 21f1ee410cSBarry Smith VecScatter *localization; /* mapping from overlapping to non-overlapping subregion; used for restrict version of algorithms */ 223c4dc4c4SMatthew Knepley VecScatter *prolongation; /* mapping from subregion to global */ 2353888de8SMatthew Knepley Vec *x,*y,*y_local; /* work vectors */ 242b691e39SMatthew Knepley IS *is; /* index set that defines each overlapping subdomain */ 252b691e39SMatthew Knepley IS *is_local; /* index set that defines each non-overlapping subdomain, may be NULL */ 264b9ad928SBarry Smith Mat *mat,*pmat; /* mat is not currently used */ 274b9ad928SBarry Smith PCASMType type; /* use reduced interpolation, restriction or both */ 28ace3abfcSBarry Smith PetscBool type_set; /* if user set this value (so won't change it for symmetric problems) */ 29ace3abfcSBarry Smith PetscBool same_local_solves; /* flag indicating whether all local solvers are same */ 30ace3abfcSBarry Smith PetscBool sort_indices; /* flag to sort subdomain indices */ 31d709ea83SDmitry Karpeev PetscBool dm_subdomains; /* whether DM is allowed to define subdomains */ 3212cd4985SMatthew G. Knepley PCCompositeType loctype; /* the type of composition for local solves */ 3380ec0b7dSPatrick Sanan MatType sub_mat_type; /* the type of Mat used for subdomain solves (can be MATSAME or NULL) */ 34fb745f2cSMatthew G. Knepley /* For multiplicative solve */ 35235cc4ceSMatthew G. Knepley Mat *lmats; /* submatrices for overlapping multiplicative (process) subdomain */ 36fb745f2cSMatthew G. Knepley Vec lx, ly; /* work vectors */ 37fb745f2cSMatthew G. Knepley IS lis; /* index set that defines each overlapping multiplicative (process) subdomain */ 385b3c0d42Seaulisa VecScatter *lprolongation; /* mapping from subregion to overlapping process subdomain */ 39347574c9Seaulisa VecScatter lrestriction; /* mapping from global to overlapping process subdomain*/ 404b9ad928SBarry Smith } PC_ASM; 414b9ad928SBarry Smith 426849ba73SBarry Smith static PetscErrorCode PCView_ASM(PC pc,PetscViewer viewer) 434b9ad928SBarry Smith { 4492bb6962SLisandro Dalcin PC_ASM *osm = (PC_ASM*)pc->data; 456849ba73SBarry Smith PetscErrorCode ierr; 4613f74950SBarry Smith PetscMPIInt rank; 474d219a6aSLisandro Dalcin PetscInt i,bsz; 48ace3abfcSBarry Smith PetscBool iascii,isstring; 494b9ad928SBarry Smith PetscViewer sviewer; 504b9ad928SBarry Smith 514b9ad928SBarry Smith PetscFunctionBegin; 52251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 53251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr); 5432077d6dSBarry Smith if (iascii) { 553d03bb48SJed Brown char overlaps[256] = "user-defined overlap",blocks[256] = "total subdomain blocks not yet set"; 568caf3d72SBarry Smith if (osm->overlap >= 0) {ierr = PetscSNPrintf(overlaps,sizeof(overlaps),"amount of overlap = %D",osm->overlap);CHKERRQ(ierr);} 578caf3d72SBarry Smith if (osm->n > 0) {ierr = PetscSNPrintf(blocks,sizeof(blocks),"total subdomain blocks = %D",osm->n);CHKERRQ(ierr);} 58efd4aadfSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %s, %s\n",blocks,overlaps);CHKERRQ(ierr); 59efd4aadfSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," restriction/interpolation type - %s\n",PCASMTypes[osm->type]);CHKERRQ(ierr); 60e64f0791SPatrick Sanan if (osm->dm_subdomains) {ierr = PetscViewerASCIIPrintf(viewer," Additive Schwarz: using DM to define subdomains\n");CHKERRQ(ierr);} 6112cd4985SMatthew G. Knepley if (osm->loctype != PC_COMPOSITE_ADDITIVE) {ierr = PetscViewerASCIIPrintf(viewer," Additive Schwarz: local solve composition type - %s\n",PCCompositeTypes[osm->loctype]);CHKERRQ(ierr);} 62ce94432eSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);CHKERRQ(ierr); 6392bb6962SLisandro Dalcin if (osm->same_local_solves) { 644d219a6aSLisandro Dalcin if (osm->ksp) { 654b9ad928SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Local solve is same for all blocks, in the following KSP and PC objects:\n");CHKERRQ(ierr); 66c5e4d11fSDmitry Karpeev ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 674d219a6aSLisandro Dalcin if (!rank) { 684b9ad928SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 6992bb6962SLisandro Dalcin ierr = KSPView(osm->ksp[0],sviewer);CHKERRQ(ierr); 704b9ad928SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 714b9ad928SBarry Smith } 72c5e4d11fSDmitry Karpeev ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 734d219a6aSLisandro Dalcin } 744b9ad928SBarry Smith } else { 75c5e4d11fSDmitry Karpeev ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 764d219a6aSLisandro Dalcin ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] number of local blocks = %D\n",(int)rank,osm->n_local_true);CHKERRQ(ierr); 774d219a6aSLisandro Dalcin ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 784b9ad928SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Local solve info for each block is in the following KSP and PC objects:\n");CHKERRQ(ierr); 794b9ad928SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 804d219a6aSLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr); 81c5e4d11fSDmitry Karpeev ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 82dd2fa690SBarry Smith for (i=0; i<osm->n_local_true; i++) { 834d219a6aSLisandro Dalcin ierr = ISGetLocalSize(osm->is[i],&bsz);CHKERRQ(ierr); 844d219a6aSLisandro Dalcin ierr = PetscViewerASCIISynchronizedPrintf(sviewer,"[%d] local block number %D, size = %D\n",(int)rank,i,bsz);CHKERRQ(ierr); 8592bb6962SLisandro Dalcin ierr = KSPView(osm->ksp[i],sviewer);CHKERRQ(ierr); 864d219a6aSLisandro Dalcin ierr = PetscViewerASCIISynchronizedPrintf(sviewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr); 874b9ad928SBarry Smith } 88c5e4d11fSDmitry Karpeev ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 894b9ad928SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 904b9ad928SBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 91c5e4d11fSDmitry Karpeev ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 924b9ad928SBarry Smith } 934b9ad928SBarry Smith } else if (isstring) { 944d219a6aSLisandro Dalcin ierr = PetscViewerStringSPrintf(viewer," blocks=%D, overlap=%D, type=%s",osm->n,osm->overlap,PCASMTypes[osm->type]);CHKERRQ(ierr); 95c5e4d11fSDmitry Karpeev ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 9692bb6962SLisandro Dalcin if (osm->ksp) {ierr = KSPView(osm->ksp[0],sviewer);CHKERRQ(ierr);} 97c5e4d11fSDmitry Karpeev ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 984b9ad928SBarry Smith } 994b9ad928SBarry Smith PetscFunctionReturn(0); 1004b9ad928SBarry Smith } 1014b9ad928SBarry Smith 10292bb6962SLisandro Dalcin static PetscErrorCode PCASMPrintSubdomains(PC pc) 10392bb6962SLisandro Dalcin { 10492bb6962SLisandro Dalcin PC_ASM *osm = (PC_ASM*)pc->data; 10592bb6962SLisandro Dalcin const char *prefix; 10692bb6962SLisandro Dalcin char fname[PETSC_MAX_PATH_LEN+1]; 107643535aeSDmitry Karpeev PetscViewer viewer, sviewer; 108643535aeSDmitry Karpeev char *s; 10992bb6962SLisandro Dalcin PetscInt i,j,nidx; 11092bb6962SLisandro Dalcin const PetscInt *idx; 111643535aeSDmitry Karpeev PetscMPIInt rank, size; 11292bb6962SLisandro Dalcin PetscErrorCode ierr; 113846783a0SBarry Smith 11492bb6962SLisandro Dalcin PetscFunctionBegin; 115ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size);CHKERRQ(ierr); 116ce94432eSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank);CHKERRQ(ierr); 11792bb6962SLisandro Dalcin ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 118c5929fdfSBarry Smith ierr = PetscOptionsGetString(NULL,prefix,"-pc_asm_print_subdomains",fname,PETSC_MAX_PATH_LEN,NULL);CHKERRQ(ierr); 11992bb6962SLisandro Dalcin if (fname[0] == 0) { ierr = PetscStrcpy(fname,"stdout");CHKERRQ(ierr); }; 120ce94432eSBarry Smith ierr = PetscViewerASCIIOpen(PetscObjectComm((PetscObject)pc),fname,&viewer);CHKERRQ(ierr); 121643535aeSDmitry Karpeev for (i=0; i<osm->n_local; i++) { 122643535aeSDmitry Karpeev if (i < osm->n_local_true) { 12392bb6962SLisandro Dalcin ierr = ISGetLocalSize(osm->is[i],&nidx);CHKERRQ(ierr); 12492bb6962SLisandro Dalcin ierr = ISGetIndices(osm->is[i],&idx);CHKERRQ(ierr); 125643535aeSDmitry Karpeev /* Print to a string viewer; no more than 15 characters per index plus 512 char for the header.*/ 126854ce69bSBarry Smith ierr = PetscMalloc1(16*(nidx+1)+512, &s);CHKERRQ(ierr); 127643535aeSDmitry Karpeev ierr = PetscViewerStringOpen(PETSC_COMM_SELF, s, 16*(nidx+1)+512, &sviewer);CHKERRQ(ierr); 128643535aeSDmitry Karpeev ierr = PetscViewerStringSPrintf(sviewer, "[%D:%D] Subdomain %D with overlap:\n", rank, size, i);CHKERRQ(ierr); 12992bb6962SLisandro Dalcin for (j=0; j<nidx; j++) { 130643535aeSDmitry Karpeev ierr = PetscViewerStringSPrintf(sviewer,"%D ",idx[j]);CHKERRQ(ierr); 13192bb6962SLisandro Dalcin } 13292bb6962SLisandro Dalcin ierr = ISRestoreIndices(osm->is[i],&idx);CHKERRQ(ierr); 133643535aeSDmitry Karpeev ierr = PetscViewerStringSPrintf(sviewer,"\n");CHKERRQ(ierr); 134643535aeSDmitry Karpeev ierr = PetscViewerDestroy(&sviewer);CHKERRQ(ierr); 135c5e4d11fSDmitry Karpeev ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 136643535aeSDmitry Karpeev ierr = PetscViewerASCIISynchronizedPrintf(viewer, s);CHKERRQ(ierr); 137643535aeSDmitry Karpeev ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 138c5e4d11fSDmitry Karpeev ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 139643535aeSDmitry Karpeev ierr = PetscFree(s);CHKERRQ(ierr); 1402b691e39SMatthew Knepley if (osm->is_local) { 141643535aeSDmitry Karpeev /* Print to a string viewer; no more than 15 characters per index plus 512 char for the header.*/ 142854ce69bSBarry Smith ierr = PetscMalloc1(16*(nidx+1)+512, &s);CHKERRQ(ierr); 143643535aeSDmitry Karpeev ierr = PetscViewerStringOpen(PETSC_COMM_SELF, s, 16*(nidx+1)+512, &sviewer);CHKERRQ(ierr); 14409d011a0SDmitry Karpeev ierr = PetscViewerStringSPrintf(sviewer, "[%D:%D] Subdomain %D without overlap:\n", rank, size, i);CHKERRQ(ierr); 1452b691e39SMatthew Knepley ierr = ISGetLocalSize(osm->is_local[i],&nidx);CHKERRQ(ierr); 1462b691e39SMatthew Knepley ierr = ISGetIndices(osm->is_local[i],&idx);CHKERRQ(ierr); 1472b691e39SMatthew Knepley for (j=0; j<nidx; j++) { 14809d011a0SDmitry Karpeev ierr = PetscViewerStringSPrintf(sviewer,"%D ",idx[j]);CHKERRQ(ierr); 1492b691e39SMatthew Knepley } 1502b691e39SMatthew Knepley ierr = ISRestoreIndices(osm->is_local[i],&idx);CHKERRQ(ierr); 15109d011a0SDmitry Karpeev ierr = PetscViewerStringSPrintf(sviewer,"\n");CHKERRQ(ierr); 152643535aeSDmitry Karpeev ierr = PetscViewerDestroy(&sviewer);CHKERRQ(ierr); 153c5e4d11fSDmitry Karpeev ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 154643535aeSDmitry Karpeev ierr = PetscViewerASCIISynchronizedPrintf(viewer, s);CHKERRQ(ierr); 155643535aeSDmitry Karpeev ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 156c5e4d11fSDmitry Karpeev ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 157643535aeSDmitry Karpeev ierr = PetscFree(s);CHKERRQ(ierr); 158643535aeSDmitry Karpeev } 1592fa5cd67SKarl Rupp } else { 160643535aeSDmitry Karpeev /* Participate in collective viewer calls. */ 161c5e4d11fSDmitry Karpeev ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 162643535aeSDmitry Karpeev ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 163c5e4d11fSDmitry Karpeev ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 164643535aeSDmitry Karpeev /* Assume either all ranks have is_local or none do. */ 165643535aeSDmitry Karpeev if (osm->is_local) { 166c5e4d11fSDmitry Karpeev ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 167643535aeSDmitry Karpeev ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 168c5e4d11fSDmitry Karpeev ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 169643535aeSDmitry Karpeev } 170fdc48646SDmitry Karpeev } 17192bb6962SLisandro Dalcin } 17292bb6962SLisandro Dalcin ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 173fcfd50ebSBarry Smith ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 17492bb6962SLisandro Dalcin PetscFunctionReturn(0); 17592bb6962SLisandro Dalcin } 17692bb6962SLisandro Dalcin 1776849ba73SBarry Smith static PetscErrorCode PCSetUp_ASM(PC pc) 1784b9ad928SBarry Smith { 1794b9ad928SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 1806849ba73SBarry Smith PetscErrorCode ierr; 181ace3abfcSBarry Smith PetscBool symset,flg; 18287e86712SBarry Smith PetscInt i,m,m_local; 1834b9ad928SBarry Smith MatReuse scall = MAT_REUSE_MATRIX; 1844b9ad928SBarry Smith IS isl; 1854b9ad928SBarry Smith KSP ksp; 1864b9ad928SBarry Smith PC subpc; 1872dcb1b2aSMatthew Knepley const char *prefix,*pprefix; 18823ce1328SBarry Smith Vec vec; 1890298fd71SBarry Smith DM *domain_dm = NULL; 1904b9ad928SBarry Smith 1914b9ad928SBarry Smith PetscFunctionBegin; 1924b9ad928SBarry Smith if (!pc->setupcalled) { 19392bb6962SLisandro Dalcin 19492bb6962SLisandro Dalcin if (!osm->type_set) { 19592bb6962SLisandro Dalcin ierr = MatIsSymmetricKnown(pc->pmat,&symset,&flg);CHKERRQ(ierr); 1962fa5cd67SKarl Rupp if (symset && flg) osm->type = PC_ASM_BASIC; 19792bb6962SLisandro Dalcin } 19892bb6962SLisandro Dalcin 1992b837212SDmitry Karpeev /* Note: if subdomains have been set either via PCASMSetTotalSubdomains() or via PCASMSetLocalSubdomains(), osm->n_local_true will not be PETSC_DECIDE */ 2002b837212SDmitry Karpeev if (osm->n_local_true == PETSC_DECIDE) { 20169ca1f37SDmitry Karpeev /* no subdomains given */ 20265db9045SDmitry Karpeev /* try pc->dm first, if allowed */ 203d709ea83SDmitry Karpeev if (osm->dm_subdomains && pc->dm) { 204feb221f8SDmitry Karpeev PetscInt num_domains, d; 205feb221f8SDmitry Karpeev char **domain_names; 2068d4ac253SDmitry Karpeev IS *inner_domain_is, *outer_domain_is; 2078d4ac253SDmitry Karpeev ierr = DMCreateDomainDecomposition(pc->dm, &num_domains, &domain_names, &inner_domain_is, &outer_domain_is, &domain_dm);CHKERRQ(ierr); 208704f0395SPatrick Sanan osm->overlap = -1; /* We do not want to increase the overlap of the IS. 209704f0395SPatrick Sanan A future improvement of this code might allow one to use 210704f0395SPatrick Sanan DM-defined subdomains and also increase the overlap, 211704f0395SPatrick Sanan but that is not currently supported */ 21269ca1f37SDmitry Karpeev if (num_domains) { 2138d4ac253SDmitry Karpeev ierr = PCASMSetLocalSubdomains(pc, num_domains, outer_domain_is, inner_domain_is);CHKERRQ(ierr); 21469ca1f37SDmitry Karpeev } 215feb221f8SDmitry Karpeev for (d = 0; d < num_domains; ++d) { 216a35b7b57SDmitry Karpeev if (domain_names) {ierr = PetscFree(domain_names[d]);CHKERRQ(ierr);} 217a35b7b57SDmitry Karpeev if (inner_domain_is) {ierr = ISDestroy(&inner_domain_is[d]);CHKERRQ(ierr);} 218a35b7b57SDmitry Karpeev if (outer_domain_is) {ierr = ISDestroy(&outer_domain_is[d]);CHKERRQ(ierr);} 219feb221f8SDmitry Karpeev } 220feb221f8SDmitry Karpeev ierr = PetscFree(domain_names);CHKERRQ(ierr); 2218d4ac253SDmitry Karpeev ierr = PetscFree(inner_domain_is);CHKERRQ(ierr); 2228d4ac253SDmitry Karpeev ierr = PetscFree(outer_domain_is);CHKERRQ(ierr); 223feb221f8SDmitry Karpeev } 2242b837212SDmitry Karpeev if (osm->n_local_true == PETSC_DECIDE) { 22569ca1f37SDmitry Karpeev /* still no subdomains; use one subdomain per processor */ 2262b837212SDmitry Karpeev osm->n_local_true = 1; 227feb221f8SDmitry Karpeev } 2282b837212SDmitry Karpeev } 2292b837212SDmitry Karpeev { /* determine the global and max number of subdomains */ 2306ac3741eSJed Brown struct {PetscInt max,sum;} inwork,outwork; 231c10200c1SHong Zhang PetscMPIInt size; 232c10200c1SHong Zhang 2336ac3741eSJed Brown inwork.max = osm->n_local_true; 2346ac3741eSJed Brown inwork.sum = osm->n_local_true; 235367daffbSBarry Smith ierr = MPIU_Allreduce(&inwork,&outwork,1,MPIU_2INT,MPIU_MAXSUM_OP,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr); 2366ac3741eSJed Brown osm->n_local = outwork.max; 2376ac3741eSJed Brown osm->n = outwork.sum; 238c10200c1SHong Zhang 239c10200c1SHong Zhang ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr); 240c10200c1SHong Zhang if (outwork.max == 1 && outwork.sum == size) { 2417dae84e0SHong Zhang /* osm->n_local_true = 1 on all processes, set this option may enable use of optimized MatCreateSubMatrices() implementation */ 242c10200c1SHong Zhang ierr = MatSetOption(pc->pmat,MAT_SUBMAT_SINGLEIS,PETSC_TRUE);CHKERRQ(ierr); 243c10200c1SHong Zhang } 2444b9ad928SBarry Smith } 24578904715SLisandro Dalcin if (!osm->is) { /* create the index sets */ 24678904715SLisandro Dalcin ierr = PCASMCreateSubdomains(pc->pmat,osm->n_local_true,&osm->is);CHKERRQ(ierr); 2474b9ad928SBarry Smith } 248f5234e35SJed Brown if (osm->n_local_true > 1 && !osm->is_local) { 249785e854fSJed Brown ierr = PetscMalloc1(osm->n_local_true,&osm->is_local);CHKERRQ(ierr); 250f5234e35SJed Brown for (i=0; i<osm->n_local_true; i++) { 251f5234e35SJed Brown if (osm->overlap > 0) { /* With positive overlap, osm->is[i] will be modified */ 252f5234e35SJed Brown ierr = ISDuplicate(osm->is[i],&osm->is_local[i]);CHKERRQ(ierr); 253f5234e35SJed Brown ierr = ISCopy(osm->is[i],osm->is_local[i]);CHKERRQ(ierr); 254f5234e35SJed Brown } else { 255f5234e35SJed Brown ierr = PetscObjectReference((PetscObject)osm->is[i]);CHKERRQ(ierr); 256f5234e35SJed Brown osm->is_local[i] = osm->is[i]; 257f5234e35SJed Brown } 258f5234e35SJed Brown } 259f5234e35SJed Brown } 26092bb6962SLisandro Dalcin ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 26190d69ab7SBarry Smith flg = PETSC_FALSE; 262c5929fdfSBarry Smith ierr = PetscOptionsGetBool(NULL,prefix,"-pc_asm_print_subdomains",&flg,NULL);CHKERRQ(ierr); 26392bb6962SLisandro Dalcin if (flg) { ierr = PCASMPrintSubdomains(pc);CHKERRQ(ierr); } 2644b9ad928SBarry Smith 2653d03bb48SJed Brown if (osm->overlap > 0) { 2664b9ad928SBarry Smith /* Extend the "overlapping" regions by a number of steps */ 26792bb6962SLisandro Dalcin ierr = MatIncreaseOverlap(pc->pmat,osm->n_local_true,osm->is,osm->overlap);CHKERRQ(ierr); 2683d03bb48SJed Brown } 2696ed231c7SMatthew Knepley if (osm->sort_indices) { 27092bb6962SLisandro Dalcin for (i=0; i<osm->n_local_true; i++) { 2714b9ad928SBarry Smith ierr = ISSort(osm->is[i]);CHKERRQ(ierr); 2722b691e39SMatthew Knepley if (osm->is_local) { 2732b691e39SMatthew Knepley ierr = ISSort(osm->is_local[i]);CHKERRQ(ierr); 2742b691e39SMatthew Knepley } 2754b9ad928SBarry Smith } 2766ed231c7SMatthew Knepley } 2774b9ad928SBarry Smith 278f6991133SBarry Smith if (!osm->ksp) { 27978904715SLisandro Dalcin /* Create the local solvers */ 280785e854fSJed Brown ierr = PetscMalloc1(osm->n_local_true,&osm->ksp);CHKERRQ(ierr); 281feb221f8SDmitry Karpeev if (domain_dm) { 282feb221f8SDmitry Karpeev ierr = PetscInfo(pc,"Setting up ASM subproblems using the embedded DM\n");CHKERRQ(ierr); 283feb221f8SDmitry Karpeev } 28492bb6962SLisandro Dalcin for (i=0; i<osm->n_local_true; i++) { 2854b9ad928SBarry Smith ierr = KSPCreate(PETSC_COMM_SELF,&ksp);CHKERRQ(ierr); 286422a814eSBarry Smith ierr = KSPSetErrorIfNotConverged(ksp,pc->erroriffailure);CHKERRQ(ierr); 2873bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ksp);CHKERRQ(ierr); 28892bb6962SLisandro Dalcin ierr = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,1);CHKERRQ(ierr); 2894b9ad928SBarry Smith ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr); 2904b9ad928SBarry Smith ierr = KSPGetPC(ksp,&subpc);CHKERRQ(ierr); 2914b9ad928SBarry Smith ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 2924b9ad928SBarry Smith ierr = KSPSetOptionsPrefix(ksp,prefix);CHKERRQ(ierr); 2934b9ad928SBarry Smith ierr = KSPAppendOptionsPrefix(ksp,"sub_");CHKERRQ(ierr); 294feb221f8SDmitry Karpeev if (domain_dm) { 295feb221f8SDmitry Karpeev ierr = KSPSetDM(ksp, domain_dm[i]);CHKERRQ(ierr); 296feb221f8SDmitry Karpeev ierr = KSPSetDMActive(ksp, PETSC_FALSE);CHKERRQ(ierr); 297feb221f8SDmitry Karpeev ierr = DMDestroy(&domain_dm[i]);CHKERRQ(ierr); 298feb221f8SDmitry Karpeev } 2994b9ad928SBarry Smith osm->ksp[i] = ksp; 3004b9ad928SBarry Smith } 301feb221f8SDmitry Karpeev if (domain_dm) { 302feb221f8SDmitry Karpeev ierr = PetscFree(domain_dm);CHKERRQ(ierr); 303feb221f8SDmitry Karpeev } 304f6991133SBarry Smith } 305347574c9Seaulisa ierr = ISConcatenate(PETSC_COMM_SELF, osm->n_local_true, osm->is, &osm->lis);CHKERRQ(ierr); 306347574c9Seaulisa ierr = ISSortRemoveDups(osm->lis);CHKERRQ(ierr); 307347574c9Seaulisa //if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE) { 308fb745f2cSMatthew G. Knepley PetscInt m; 309fb745f2cSMatthew G. Knepley ierr = ISConcatenate(PETSC_COMM_SELF, osm->n_local_true, osm->is, &osm->lis);CHKERRQ(ierr); 310fb745f2cSMatthew G. Knepley ierr = ISSortRemoveDups(osm->lis);CHKERRQ(ierr); 311fb745f2cSMatthew G. Knepley ierr = ISGetLocalSize(osm->lis, &m);CHKERRQ(ierr); 312fb745f2cSMatthew G. Knepley ierr = VecCreateSeq(PETSC_COMM_SELF, m, &osm->lx);CHKERRQ(ierr); 313fb745f2cSMatthew G. Knepley ierr = VecDuplicate(osm->lx, &osm->ly);CHKERRQ(ierr); 314347574c9Seaulisa //} 3154b9ad928SBarry Smith scall = MAT_INITIAL_MATRIX; 3164b9ad928SBarry Smith } else { 3174b9ad928SBarry Smith /* 3184b9ad928SBarry Smith Destroy the blocks from the previous iteration 3194b9ad928SBarry Smith */ 320e09e08ccSBarry Smith if (pc->flag == DIFFERENT_NONZERO_PATTERN) { 3214b9ad928SBarry Smith ierr = MatDestroyMatrices(osm->n_local_true,&osm->pmat);CHKERRQ(ierr); 3224b9ad928SBarry Smith scall = MAT_INITIAL_MATRIX; 3234b9ad928SBarry Smith } 3244b9ad928SBarry Smith } 3254b9ad928SBarry Smith 32692bb6962SLisandro Dalcin /* 32792bb6962SLisandro Dalcin Extract out the submatrices 32892bb6962SLisandro Dalcin */ 3297dae84e0SHong Zhang ierr = MatCreateSubMatrices(pc->pmat,osm->n_local_true,osm->is,osm->is,scall,&osm->pmat);CHKERRQ(ierr); 33092bb6962SLisandro Dalcin if (scall == MAT_INITIAL_MATRIX) { 33192bb6962SLisandro Dalcin ierr = PetscObjectGetOptionsPrefix((PetscObject)pc->pmat,&pprefix);CHKERRQ(ierr); 33292bb6962SLisandro Dalcin for (i=0; i<osm->n_local_true; i++) { 3333bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)osm->pmat[i]);CHKERRQ(ierr); 33478904715SLisandro Dalcin ierr = PetscObjectSetOptionsPrefix((PetscObject)osm->pmat[i],pprefix);CHKERRQ(ierr); 33592bb6962SLisandro Dalcin } 33692bb6962SLisandro Dalcin } 33780ec0b7dSPatrick Sanan 33880ec0b7dSPatrick Sanan /* Convert the types of the submatrices (if needbe) */ 33980ec0b7dSPatrick Sanan if (osm->sub_mat_type) { 34080ec0b7dSPatrick Sanan for (i=0; i<osm->n_local_true; i++) { 34180ec0b7dSPatrick Sanan ierr = MatConvert(osm->pmat[i],osm->sub_mat_type,MAT_INPLACE_MATRIX,&(osm->pmat[i]));CHKERRQ(ierr); 34280ec0b7dSPatrick Sanan } 34380ec0b7dSPatrick Sanan } 34480ec0b7dSPatrick Sanan 34580ec0b7dSPatrick Sanan if(!pc->setupcalled){ 34680ec0b7dSPatrick Sanan /* Create the local work vectors (from the local matrices) and scatter contexts */ 34780ec0b7dSPatrick Sanan ierr = MatCreateVecs(pc->pmat,&vec,0);CHKERRQ(ierr); 34880ec0b7dSPatrick Sanan ierr = PetscMalloc1(osm->n_local,&osm->restriction);CHKERRQ(ierr); 349f1ee410cSBarry Smith if (osm->is_local && (osm->type == PC_ASM_INTERPOLATE || osm->type == PC_ASM_NONE )) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Cannot use interpolate or none PCASMType if is_local was provided to PCASMSetLocalSubdomains()"); 350*d8b3b5e3Seaulisa if (osm->is_local && osm->type == PC_ASM_RESTRICT && osm->loctype == PC_COMPOSITE_ADDITIVE) {ierr = PetscMalloc1(osm->n_local,&osm->localization);CHKERRQ(ierr);} 35180ec0b7dSPatrick Sanan ierr = PetscMalloc1(osm->n_local,&osm->prolongation);CHKERRQ(ierr); 352347574c9Seaulisa //if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE) 3535b3c0d42Seaulisa ierr = PetscMalloc1(osm->n_local,&osm->lprolongation);CHKERRQ(ierr); 3545b3c0d42Seaulisa 35580ec0b7dSPatrick Sanan ierr = PetscMalloc1(osm->n_local,&osm->x);CHKERRQ(ierr); 35680ec0b7dSPatrick Sanan ierr = PetscMalloc1(osm->n_local,&osm->y);CHKERRQ(ierr); 35780ec0b7dSPatrick Sanan ierr = PetscMalloc1(osm->n_local,&osm->y_local);CHKERRQ(ierr); 358347574c9Seaulisa 359347574c9Seaulisa ierr = ISGetLocalSize(osm->lis,&m);CHKERRQ(ierr); 360347574c9Seaulisa ierr = ISCreateStride(PETSC_COMM_SELF,m,0,1,&isl);CHKERRQ(ierr); 361347574c9Seaulisa ierr = VecScatterCreate(vec,osm->lis,osm->lx,isl,&osm->lrestriction);CHKERRQ(ierr); 362347574c9Seaulisa ierr = ISDestroy(&isl);CHKERRQ(ierr); 363347574c9Seaulisa 364347574c9Seaulisa 36580ec0b7dSPatrick Sanan for (i=0; i<osm->n_local_true; ++i) { 36680ec0b7dSPatrick Sanan ierr = ISGetLocalSize(osm->is[i],&m);CHKERRQ(ierr); 36780ec0b7dSPatrick Sanan ierr = MatCreateVecs(osm->pmat[i],&osm->x[i],NULL);CHKERRQ(ierr); 36880ec0b7dSPatrick Sanan ierr = ISCreateStride(PETSC_COMM_SELF,m,0,1,&isl);CHKERRQ(ierr); 36980ec0b7dSPatrick Sanan ierr = VecScatterCreate(vec,osm->is[i],osm->x[i],isl,&osm->restriction[i]);CHKERRQ(ierr); 37080ec0b7dSPatrick Sanan ierr = ISDestroy(&isl);CHKERRQ(ierr); 37180ec0b7dSPatrick Sanan ierr = VecDuplicate(osm->x[i],&osm->y[i]);CHKERRQ(ierr); 372f1ee410cSBarry Smith if (osm->localization) { 37380ec0b7dSPatrick Sanan ISLocalToGlobalMapping ltog; 37480ec0b7dSPatrick Sanan IS isll; 37580ec0b7dSPatrick Sanan const PetscInt *idx_local; 37680ec0b7dSPatrick Sanan PetscInt *idx,nout; 37780ec0b7dSPatrick Sanan 37880ec0b7dSPatrick Sanan ierr = ISLocalToGlobalMappingCreateIS(osm->is[i],<og);CHKERRQ(ierr); 37980ec0b7dSPatrick Sanan ierr = ISGetLocalSize(osm->is_local[i],&m_local);CHKERRQ(ierr); 38080ec0b7dSPatrick Sanan ierr = ISGetIndices(osm->is_local[i], &idx_local);CHKERRQ(ierr); 38180ec0b7dSPatrick Sanan ierr = PetscMalloc1(m_local,&idx);CHKERRQ(ierr); 38280ec0b7dSPatrick Sanan ierr = ISGlobalToLocalMappingApply(ltog,IS_GTOLM_DROP,m_local,idx_local,&nout,idx);CHKERRQ(ierr); 38380ec0b7dSPatrick Sanan ierr = ISLocalToGlobalMappingDestroy(<og);CHKERRQ(ierr); 38480ec0b7dSPatrick Sanan if (nout != m_local) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"is_local not a subset of is"); 38580ec0b7dSPatrick Sanan ierr = ISRestoreIndices(osm->is_local[i], &idx_local);CHKERRQ(ierr); 38680ec0b7dSPatrick Sanan ierr = ISCreateGeneral(PETSC_COMM_SELF,m_local,idx,PETSC_OWN_POINTER,&isll);CHKERRQ(ierr); 38780ec0b7dSPatrick Sanan ierr = ISCreateStride(PETSC_COMM_SELF,m_local,0,1,&isl);CHKERRQ(ierr); 38880ec0b7dSPatrick Sanan ierr = VecCreateSeq(PETSC_COMM_SELF,m_local,&osm->y_local[i]);CHKERRQ(ierr); 389*d8b3b5e3Seaulisa //ierr = VecScatterCreate(osm->y[i],isll,osm->y_local[i],isl,&osm->localization[i]);CHKERRQ(ierr); 39080ec0b7dSPatrick Sanan ierr = ISDestroy(&isll);CHKERRQ(ierr); 39180ec0b7dSPatrick Sanan 39280ec0b7dSPatrick Sanan ierr = VecScatterCreate(vec,osm->is_local[i],osm->y_local[i],isl,&osm->prolongation[i]);CHKERRQ(ierr); 39380ec0b7dSPatrick Sanan ierr = ISDestroy(&isl);CHKERRQ(ierr); 39480ec0b7dSPatrick Sanan } else { 39580ec0b7dSPatrick Sanan osm->y_local[i] = osm->y[i]; 39680ec0b7dSPatrick Sanan ierr = PetscObjectReference((PetscObject) osm->y[i]);CHKERRQ(ierr); 39780ec0b7dSPatrick Sanan osm->prolongation[i] = osm->restriction[i]; 39880ec0b7dSPatrick Sanan ierr = PetscObjectReference((PetscObject) osm->restriction[i]);CHKERRQ(ierr); 39980ec0b7dSPatrick Sanan } 4005b3c0d42Seaulisa 401*d8b3b5e3Seaulisa 4025b3c0d42Seaulisa ISLocalToGlobalMapping ltog; 4035b3c0d42Seaulisa IS isll; 4045b3c0d42Seaulisa const PetscInt *idx_is; 4055b3c0d42Seaulisa PetscInt *idx_lis,nout; 4065b3c0d42Seaulisa 4075b3c0d42Seaulisa ierr = ISLocalToGlobalMappingCreateIS(osm->lis,<og);CHKERRQ(ierr); 4085b3c0d42Seaulisa ierr = ISGetLocalSize(osm->is[i],&m);CHKERRQ(ierr); 4095b3c0d42Seaulisa ierr = ISGetIndices(osm->is[i], &idx_is);CHKERRQ(ierr); 4105b3c0d42Seaulisa ierr = PetscMalloc1(m,&idx_lis);CHKERRQ(ierr); 4115b3c0d42Seaulisa ierr = ISGlobalToLocalMappingApply(ltog,IS_GTOLM_DROP,m,idx_is,&nout,idx_lis);CHKERRQ(ierr); 4125b3c0d42Seaulisa if (nout != m) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"is not a subset of lis"); 4135b3c0d42Seaulisa ierr = ISRestoreIndices(osm->is[i], &idx_is);CHKERRQ(ierr); 4145b3c0d42Seaulisa ierr = ISCreateGeneral(PETSC_COMM_SELF,m,idx_lis,PETSC_OWN_POINTER,&isll);CHKERRQ(ierr); 415*d8b3b5e3Seaulisa ierr = ISLocalToGlobalMappingDestroy(<og);CHKERRQ(ierr); 4165b3c0d42Seaulisa ierr = ISCreateStride(PETSC_COMM_SELF,m,0,1,&isl);CHKERRQ(ierr); 4175b3c0d42Seaulisa ierr = VecScatterCreate(osm->ly,isll,osm->y[i],isl,&osm->lprolongation[i]);CHKERRQ(ierr); 4185b3c0d42Seaulisa ierr = ISDestroy(&isll);CHKERRQ(ierr); 4195b3c0d42Seaulisa ierr = ISDestroy(&isl);CHKERRQ(ierr); 420*d8b3b5e3Seaulisa if (osm->localization) { 421*d8b3b5e3Seaulisa ISLocalToGlobalMapping ltog; 422*d8b3b5e3Seaulisa IS isll,isll_local; 423*d8b3b5e3Seaulisa const PetscInt *idx_local; 424*d8b3b5e3Seaulisa PetscInt *idx1, *idx2, nout; 425*d8b3b5e3Seaulisa 426*d8b3b5e3Seaulisa ierr = ISGetLocalSize(osm->is_local[i],&m_local);CHKERRQ(ierr); 427*d8b3b5e3Seaulisa ierr = ISGetIndices(osm->is_local[i], &idx_local);CHKERRQ(ierr); 428*d8b3b5e3Seaulisa 429*d8b3b5e3Seaulisa ierr = ISLocalToGlobalMappingCreateIS(osm->is[i],<og);CHKERRQ(ierr); 430*d8b3b5e3Seaulisa ierr = PetscMalloc1(m_local,&idx1);CHKERRQ(ierr); 431*d8b3b5e3Seaulisa ierr = ISGlobalToLocalMappingApply(ltog,IS_GTOLM_DROP,m_local,idx_local,&nout,idx1);CHKERRQ(ierr); 432*d8b3b5e3Seaulisa ierr = ISLocalToGlobalMappingDestroy(<og);CHKERRQ(ierr); 433*d8b3b5e3Seaulisa if (nout != m_local) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"is_local not a subset of is"); 434*d8b3b5e3Seaulisa ierr = ISCreateGeneral(PETSC_COMM_SELF,m_local,idx1,PETSC_OWN_POINTER,&isll);CHKERRQ(ierr); 435*d8b3b5e3Seaulisa 436*d8b3b5e3Seaulisa ierr = ISLocalToGlobalMappingCreateIS(osm->lis,<og);CHKERRQ(ierr); 437*d8b3b5e3Seaulisa ierr = PetscMalloc1(m_local,&idx2);CHKERRQ(ierr); 438*d8b3b5e3Seaulisa ierr = ISGlobalToLocalMappingApply(ltog,IS_GTOLM_DROP,m_local,idx_local,&nout,idx2);CHKERRQ(ierr); 439*d8b3b5e3Seaulisa ierr = ISLocalToGlobalMappingDestroy(<og);CHKERRQ(ierr); 440*d8b3b5e3Seaulisa if (nout != m_local) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"is_local not a subset of lis"); 441*d8b3b5e3Seaulisa ierr = ISCreateGeneral(PETSC_COMM_SELF,m_local,idx2,PETSC_OWN_POINTER,&isll_local);CHKERRQ(ierr); 442*d8b3b5e3Seaulisa 443*d8b3b5e3Seaulisa ierr = ISRestoreIndices(osm->is_local[i], &idx_local);CHKERRQ(ierr); 444*d8b3b5e3Seaulisa ierr = VecScatterCreate(osm->y[i],isll,osm->ly,isll_local,&osm->localization[i]);CHKERRQ(ierr); 445*d8b3b5e3Seaulisa 446*d8b3b5e3Seaulisa ierr = ISDestroy(&isll);CHKERRQ(ierr); 447*d8b3b5e3Seaulisa ierr = ISDestroy(&isll_local);CHKERRQ(ierr); 448*d8b3b5e3Seaulisa } 449*d8b3b5e3Seaulisa 4505b3c0d42Seaulisa 45180ec0b7dSPatrick Sanan } 45280ec0b7dSPatrick Sanan for (i=osm->n_local_true; i<osm->n_local; i++) { 45380ec0b7dSPatrick Sanan ierr = VecCreateSeq(PETSC_COMM_SELF,0,&osm->x[i]);CHKERRQ(ierr); 45480ec0b7dSPatrick Sanan ierr = VecDuplicate(osm->x[i],&osm->y[i]);CHKERRQ(ierr); 45580ec0b7dSPatrick Sanan ierr = VecDuplicate(osm->x[i],&osm->y_local[i]);CHKERRQ(ierr); 45680ec0b7dSPatrick Sanan ierr = ISCreateStride(PETSC_COMM_SELF,0,0,1,&isl);CHKERRQ(ierr); 45780ec0b7dSPatrick Sanan ierr = VecScatterCreate(vec,isl,osm->x[i],isl,&osm->restriction[i]);CHKERRQ(ierr); 458f1ee410cSBarry Smith if (osm->localization) { 45980ec0b7dSPatrick Sanan ierr = VecScatterCreate(osm->y[i],isl,osm->y_local[i],isl,&osm->localization[i]);CHKERRQ(ierr); 46080ec0b7dSPatrick Sanan ierr = VecScatterCreate(vec,isl,osm->x[i],isl,&osm->prolongation[i]);CHKERRQ(ierr); 46180ec0b7dSPatrick Sanan } else { 46280ec0b7dSPatrick Sanan osm->prolongation[i] = osm->restriction[i]; 46380ec0b7dSPatrick Sanan ierr = PetscObjectReference((PetscObject) osm->restriction[i]);CHKERRQ(ierr); 46480ec0b7dSPatrick Sanan } 46580ec0b7dSPatrick Sanan ierr = ISDestroy(&isl);CHKERRQ(ierr); 46680ec0b7dSPatrick Sanan } 46780ec0b7dSPatrick Sanan ierr = VecDestroy(&vec);CHKERRQ(ierr); 46880ec0b7dSPatrick Sanan } 46980ec0b7dSPatrick Sanan 470fb745f2cSMatthew G. Knepley if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE) { 471235cc4ceSMatthew G. Knepley IS *cis; 472235cc4ceSMatthew G. Knepley PetscInt c; 473235cc4ceSMatthew G. Knepley 474235cc4ceSMatthew G. Knepley ierr = PetscMalloc1(osm->n_local_true, &cis);CHKERRQ(ierr); 475235cc4ceSMatthew G. Knepley for (c = 0; c < osm->n_local_true; ++c) cis[c] = osm->lis; 4767dae84e0SHong Zhang ierr = MatCreateSubMatrices(pc->pmat, osm->n_local_true, osm->is, cis, scall, &osm->lmats);CHKERRQ(ierr); 477235cc4ceSMatthew G. Knepley ierr = PetscFree(cis);CHKERRQ(ierr); 478fb745f2cSMatthew G. Knepley } 4794b9ad928SBarry Smith 4804b9ad928SBarry Smith /* Return control to the user so that the submatrices can be modified (e.g., to apply 4814b9ad928SBarry Smith different boundary conditions for the submatrices than for the global problem) */ 48292bb6962SLisandro Dalcin ierr = PCModifySubMatrices(pc,osm->n_local_true,osm->is,osm->is,osm->pmat,pc->modifysubmatricesP);CHKERRQ(ierr); 4834b9ad928SBarry Smith 48492bb6962SLisandro Dalcin /* 48592bb6962SLisandro Dalcin Loop over subdomains putting them into local ksp 48692bb6962SLisandro Dalcin */ 48792bb6962SLisandro Dalcin for (i=0; i<osm->n_local_true; i++) { 48823ee1639SBarry Smith ierr = KSPSetOperators(osm->ksp[i],osm->pmat[i],osm->pmat[i]);CHKERRQ(ierr); 48992bb6962SLisandro Dalcin if (!pc->setupcalled) { 490bf108f30SBarry Smith ierr = KSPSetFromOptions(osm->ksp[i]);CHKERRQ(ierr); 4914b9ad928SBarry Smith } 49292bb6962SLisandro Dalcin } 4934b9ad928SBarry Smith PetscFunctionReturn(0); 4944b9ad928SBarry Smith } 4954b9ad928SBarry Smith 4966849ba73SBarry Smith static PetscErrorCode PCSetUpOnBlocks_ASM(PC pc) 4974b9ad928SBarry Smith { 4984b9ad928SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 4996849ba73SBarry Smith PetscErrorCode ierr; 50013f74950SBarry Smith PetscInt i; 501539c167fSBarry Smith KSPConvergedReason reason; 5024b9ad928SBarry Smith 5034b9ad928SBarry Smith PetscFunctionBegin; 5044b9ad928SBarry Smith for (i=0; i<osm->n_local_true; i++) { 5054b9ad928SBarry Smith ierr = KSPSetUp(osm->ksp[i]);CHKERRQ(ierr); 506539c167fSBarry Smith ierr = KSPGetConvergedReason(osm->ksp[i],&reason);CHKERRQ(ierr); 507539c167fSBarry Smith if (reason == KSP_DIVERGED_PCSETUP_FAILED) { 508261222b2SHong Zhang pc->failedreason = PC_SUBPC_ERROR; 509e0eafd54SHong Zhang } 5104b9ad928SBarry Smith } 5114b9ad928SBarry Smith PetscFunctionReturn(0); 5124b9ad928SBarry Smith } 5134b9ad928SBarry Smith 5146849ba73SBarry Smith static PetscErrorCode PCApply_ASM(PC pc,Vec x,Vec y) 5154b9ad928SBarry Smith { 5164b9ad928SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 5176849ba73SBarry Smith PetscErrorCode ierr; 51813f74950SBarry Smith PetscInt i,n_local = osm->n_local,n_local_true = osm->n_local_true; 5194b9ad928SBarry Smith ScatterMode forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE; 5204b9ad928SBarry Smith 5214b9ad928SBarry Smith PetscFunctionBegin; 5224b9ad928SBarry Smith /* 5234b9ad928SBarry Smith Support for limiting the restriction or interpolation to only local 5244b9ad928SBarry Smith subdomain values (leaving the other values 0). 5254b9ad928SBarry Smith */ 5264b9ad928SBarry Smith if (!(osm->type & PC_ASM_RESTRICT)) { 5274b9ad928SBarry Smith forward = SCATTER_FORWARD_LOCAL; 5284b9ad928SBarry Smith /* have to zero the work RHS since scatter may leave some slots empty */ 529f21eeb63SBarry Smith for (i=0; i<n_local_true; i++) { 5303d03bb48SJed Brown ierr = VecZeroEntries(osm->x[i]);CHKERRQ(ierr); 5314b9ad928SBarry Smith } 5324b9ad928SBarry Smith } 533347574c9Seaulisa if (!(osm->type & PC_ASM_INTERPOLATE)) { 534347574c9Seaulisa reverse = SCATTER_REVERSE_LOCAL; 535347574c9Seaulisa } 5364b9ad928SBarry Smith 537347574c9Seaulisa if(osm->loctype == PC_COMPOSITE_MULTIPLICATIVE || osm->loctype == PC_COMPOSITE_ADDITIVE){ 538*d8b3b5e3Seaulisa //zero the global and the local solutions 53912cd4985SMatthew G. Knepley ierr = VecZeroEntries(y);CHKERRQ(ierr); 5405b3c0d42Seaulisa ierr = VecSet(osm->ly, 0.0);CHKERRQ(ierr); 541347574c9Seaulisa 542*d8b3b5e3Seaulisa // Copy the global RHS to local RHS including the ghost nodes 543347574c9Seaulisa ierr = VecScatterBegin(osm->lrestriction, x, osm->lx, INSERT_VALUES, forward);CHKERRQ(ierr); 544347574c9Seaulisa ierr = VecScatterEnd(osm->lrestriction, x, osm->lx, INSERT_VALUES, forward);CHKERRQ(ierr); 545347574c9Seaulisa 546*d8b3b5e3Seaulisa //Restrict local RHS to the overlapping 0-block RHS 547347574c9Seaulisa ierr = VecScatterBegin(osm->lprolongation[0], osm->lx, osm->x[0], INSERT_VALUES, forward);CHKERRQ(ierr); 548347574c9Seaulisa ierr = VecScatterEnd(osm->lprolongation[0], osm->lx, osm->x[0], INSERT_VALUES, forward);CHKERRQ(ierr); 549*d8b3b5e3Seaulisa 55012cd4985SMatthew G. Knepley /* do the local solves */ 55112cd4985SMatthew G. Knepley for (i = 0; i < n_local_true; ++i) { 552347574c9Seaulisa 553*d8b3b5e3Seaulisa // solve the overlapping i-block 55412cd4985SMatthew G. Knepley ierr = KSPSolve(osm->ksp[i], osm->x[i], osm->y[i]);CHKERRQ(ierr); 555*d8b3b5e3Seaulisa 556*d8b3b5e3Seaulisa if (osm->localization) { //interpolate the non-overalapping i-block solution to the local solution 557*d8b3b5e3Seaulisa ierr = VecScatterBegin(osm->localization[i], osm->y[i], osm->ly, ADD_VALUES, forward);CHKERRQ(ierr); 558*d8b3b5e3Seaulisa ierr = VecScatterEnd(osm->localization[i], osm->y[i], osm->ly, ADD_VALUES, forward);CHKERRQ(ierr); 559*d8b3b5e3Seaulisa } 560*d8b3b5e3Seaulisa else{ //interpolate the overalapping i-block solution to the local solution 5615b3c0d42Seaulisa ierr = VecScatterBegin(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, reverse);CHKERRQ(ierr); 5625b3c0d42Seaulisa ierr = VecScatterEnd(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, reverse);CHKERRQ(ierr); 563*d8b3b5e3Seaulisa } 564347574c9Seaulisa 565347574c9Seaulisa if (i < n_local_true-1) { 566*d8b3b5e3Seaulisa // Restrict local RHS to the overlapping (i+1)-block RHS 567347574c9Seaulisa ierr = VecScatterBegin(osm->lprolongation[i+1], osm->lx, osm->x[i+1], INSERT_VALUES, forward);CHKERRQ(ierr); 568347574c9Seaulisa ierr = VecScatterEnd(osm->lprolongation[i+1], osm->lx, osm->x[i+1], INSERT_VALUES, forward);CHKERRQ(ierr); 569347574c9Seaulisa 570347574c9Seaulisa if ( osm->loctype == PC_COMPOSITE_MULTIPLICATIVE){ 571*d8b3b5e3Seaulisa // udpdate the overlapping (i+1)-block RHS using the current local solution 572347574c9Seaulisa ierr = MatMult(osm->lmats[i+1], osm->ly, osm->y[i+1]);CHKERRQ(ierr); 573347574c9Seaulisa ierr = VecAXPBY(osm->x[i+1],-1.,1., osm->y[i+1]); CHKERRQ(ierr); 5747c3d802fSMatthew G. Knepley } 57512cd4985SMatthew G. Knepley } 57612cd4985SMatthew G. Knepley } 577*d8b3b5e3Seaulisa // Add the local solution to the global solution including the ghost nodes 578347574c9Seaulisa ierr = VecScatterBegin(osm->lrestriction, osm->ly, y, ADD_VALUES, reverse);CHKERRQ(ierr); 579347574c9Seaulisa ierr = VecScatterEnd(osm->lrestriction, osm->ly, y, ADD_VALUES, reverse);CHKERRQ(ierr); 580347574c9Seaulisa }else{ 581347574c9Seaulisa SETERRQ1(PetscObjectComm((PetscObject) pc), PETSC_ERR_ARG_WRONG, "Invalid local composition type: %s", PCCompositeTypes[osm->loctype]); 58212cd4985SMatthew G. Knepley } 5834b9ad928SBarry Smith PetscFunctionReturn(0); 5844b9ad928SBarry Smith } 5854b9ad928SBarry Smith 5866849ba73SBarry Smith static PetscErrorCode PCApplyTranspose_ASM(PC pc,Vec x,Vec y) 5874b9ad928SBarry Smith { 5884b9ad928SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 5896849ba73SBarry Smith PetscErrorCode ierr; 59013f74950SBarry Smith PetscInt i,n_local = osm->n_local,n_local_true = osm->n_local_true; 5914b9ad928SBarry Smith ScatterMode forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE; 5924b9ad928SBarry Smith 5934b9ad928SBarry Smith PetscFunctionBegin; 5944b9ad928SBarry Smith /* 5954b9ad928SBarry Smith Support for limiting the restriction or interpolation to only local 5964b9ad928SBarry Smith subdomain values (leaving the other values 0). 5974b9ad928SBarry Smith 5984b9ad928SBarry Smith Note: these are reversed from the PCApply_ASM() because we are applying the 5994b9ad928SBarry Smith transpose of the three terms 6004b9ad928SBarry Smith */ 601*d8b3b5e3Seaulisa 6024b9ad928SBarry Smith if (!(osm->type & PC_ASM_INTERPOLATE)) { 6034b9ad928SBarry Smith forward = SCATTER_FORWARD_LOCAL; 6044b9ad928SBarry Smith /* have to zero the work RHS since scatter may leave some slots empty */ 605f21eeb63SBarry Smith for (i=0; i<n_local_true; i++) { 60610bd88b9SJed Brown ierr = VecZeroEntries(osm->x[i]);CHKERRQ(ierr); 6074b9ad928SBarry Smith } 6084b9ad928SBarry Smith } 6092fa5cd67SKarl Rupp if (!(osm->type & PC_ASM_RESTRICT)) reverse = SCATTER_REVERSE_LOCAL; 6104b9ad928SBarry Smith 611*d8b3b5e3Seaulisa //zero the global and the local solutions 61210bd88b9SJed Brown ierr = VecZeroEntries(y);CHKERRQ(ierr); 613*d8b3b5e3Seaulisa ierr = VecSet(osm->ly, 0.0);CHKERRQ(ierr); 614*d8b3b5e3Seaulisa 615*d8b3b5e3Seaulisa // Copy the global RHS to local RHS including the ghost nodes 616*d8b3b5e3Seaulisa ierr = VecScatterBegin(osm->lrestriction, x, osm->lx, INSERT_VALUES, forward);CHKERRQ(ierr); 617*d8b3b5e3Seaulisa ierr = VecScatterEnd(osm->lrestriction, x, osm->lx, INSERT_VALUES, forward);CHKERRQ(ierr); 618*d8b3b5e3Seaulisa 619*d8b3b5e3Seaulisa //Restrict local RHS to the overlapping 0-block RHS 620*d8b3b5e3Seaulisa ierr = VecScatterBegin(osm->lprolongation[0], osm->lx, osm->x[0], INSERT_VALUES, forward);CHKERRQ(ierr); 621*d8b3b5e3Seaulisa ierr = VecScatterEnd(osm->lprolongation[0], osm->lx, osm->x[0], INSERT_VALUES, forward);CHKERRQ(ierr); 622*d8b3b5e3Seaulisa 6234b9ad928SBarry Smith /* do the local solves */ 624*d8b3b5e3Seaulisa for (i = 0; i < n_local_true; ++i) { 625*d8b3b5e3Seaulisa 626*d8b3b5e3Seaulisa // solve the overlapping i-block 627*d8b3b5e3Seaulisa ierr = KSPSolve(osm->ksp[i], osm->x[i], osm->y[i]);CHKERRQ(ierr); 628*d8b3b5e3Seaulisa 629*d8b3b5e3Seaulisa if (osm->localization) { //interpolate the non-overalapping i-block solution to the local solution 630*d8b3b5e3Seaulisa ierr = VecScatterBegin(osm->localization[i], osm->y[i], osm->ly, ADD_VALUES, forward);CHKERRQ(ierr); 631*d8b3b5e3Seaulisa ierr = VecScatterEnd(osm->localization[i], osm->y[i], osm->ly, ADD_VALUES, forward);CHKERRQ(ierr); 632b9845e0eSMatthew Knepley } 633*d8b3b5e3Seaulisa else{ //interpolate the overalapping i-block solution to the local solution 634*d8b3b5e3Seaulisa ierr = VecScatterBegin(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, reverse);CHKERRQ(ierr); 635*d8b3b5e3Seaulisa ierr = VecScatterEnd(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, reverse);CHKERRQ(ierr); 6364b9ad928SBarry Smith } 637*d8b3b5e3Seaulisa 638*d8b3b5e3Seaulisa if (i < n_local_true-1) { 639*d8b3b5e3Seaulisa // Restrict local RHS to the overlapping (i+1)-block RHS 640*d8b3b5e3Seaulisa ierr = VecScatterBegin(osm->lprolongation[i+1], osm->lx, osm->x[i+1], INSERT_VALUES, forward);CHKERRQ(ierr); 641*d8b3b5e3Seaulisa ierr = VecScatterEnd(osm->lprolongation[i+1], osm->lx, osm->x[i+1], INSERT_VALUES, forward);CHKERRQ(ierr); 6424b9ad928SBarry Smith } 6434b9ad928SBarry Smith } 644*d8b3b5e3Seaulisa // Add the local solution to the global solution including the ghost nodes 645*d8b3b5e3Seaulisa ierr = VecScatterBegin(osm->lrestriction, osm->ly, y, ADD_VALUES, reverse);CHKERRQ(ierr); 646*d8b3b5e3Seaulisa ierr = VecScatterEnd(osm->lrestriction, osm->ly, y, ADD_VALUES, reverse);CHKERRQ(ierr); 647*d8b3b5e3Seaulisa 6484b9ad928SBarry Smith PetscFunctionReturn(0); 649*d8b3b5e3Seaulisa 6504b9ad928SBarry Smith } 6514b9ad928SBarry Smith 652e91c6855SBarry Smith static PetscErrorCode PCReset_ASM(PC pc) 6534b9ad928SBarry Smith { 6544b9ad928SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 6556849ba73SBarry Smith PetscErrorCode ierr; 65613f74950SBarry Smith PetscInt i; 6574b9ad928SBarry Smith 6584b9ad928SBarry Smith PetscFunctionBegin; 65992bb6962SLisandro Dalcin if (osm->ksp) { 66092bb6962SLisandro Dalcin for (i=0; i<osm->n_local_true; i++) { 661e91c6855SBarry Smith ierr = KSPReset(osm->ksp[i]);CHKERRQ(ierr); 66292bb6962SLisandro Dalcin } 66392bb6962SLisandro Dalcin } 664e09e08ccSBarry Smith if (osm->pmat) { 66592bb6962SLisandro Dalcin if (osm->n_local_true > 0) { 66630a70a9aSHong Zhang ierr = MatDestroySubMatrices(osm->n_local_true,&osm->pmat);CHKERRQ(ierr); 66792bb6962SLisandro Dalcin } 66892bb6962SLisandro Dalcin } 6692b691e39SMatthew Knepley if (osm->restriction) { 670347574c9Seaulisa ierr = VecScatterDestroy(&osm->lrestriction);CHKERRQ(ierr); 6714b9ad928SBarry Smith for (i=0; i<osm->n_local; i++) { 672fcfd50ebSBarry Smith ierr = VecScatterDestroy(&osm->restriction[i]);CHKERRQ(ierr); 673fcfd50ebSBarry Smith if (osm->localization) {ierr = VecScatterDestroy(&osm->localization[i]);CHKERRQ(ierr);} 674fcfd50ebSBarry Smith ierr = VecScatterDestroy(&osm->prolongation[i]);CHKERRQ(ierr); 675347574c9Seaulisa // if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE && i < osm->n_local_true){// - 1){ 676347574c9Seaulisa if (i < osm->n_local_true){// - 1){ 6775b3c0d42Seaulisa ierr = VecScatterDestroy(&osm->lprolongation[i]);CHKERRQ(ierr); 6785b3c0d42Seaulisa } 679347574c9Seaulisa 680fcfd50ebSBarry Smith ierr = VecDestroy(&osm->x[i]);CHKERRQ(ierr); 681fcfd50ebSBarry Smith ierr = VecDestroy(&osm->y[i]);CHKERRQ(ierr); 682fcfd50ebSBarry Smith ierr = VecDestroy(&osm->y_local[i]);CHKERRQ(ierr); 6834b9ad928SBarry Smith } 6842b691e39SMatthew Knepley ierr = PetscFree(osm->restriction);CHKERRQ(ierr); 685b9845e0eSMatthew Knepley if (osm->localization) {ierr = PetscFree(osm->localization);CHKERRQ(ierr);} 6862b691e39SMatthew Knepley ierr = PetscFree(osm->prolongation);CHKERRQ(ierr); 687347574c9Seaulisa //if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE){ 6885b3c0d42Seaulisa ierr = PetscFree(osm->lprolongation);CHKERRQ(ierr); 689347574c9Seaulisa //} 69005b42c5fSBarry Smith ierr = PetscFree(osm->x);CHKERRQ(ierr); 69178904715SLisandro Dalcin ierr = PetscFree(osm->y);CHKERRQ(ierr); 69253888de8SMatthew Knepley ierr = PetscFree(osm->y_local);CHKERRQ(ierr); 69392bb6962SLisandro Dalcin } 6942b691e39SMatthew Knepley ierr = PCASMDestroySubdomains(osm->n_local_true,osm->is,osm->is_local);CHKERRQ(ierr); 695fb745f2cSMatthew G. Knepley ierr = ISDestroy(&osm->lis);CHKERRQ(ierr); 696fb745f2cSMatthew G. Knepley ierr = VecDestroy(&osm->lx);CHKERRQ(ierr); 697fb745f2cSMatthew G. Knepley ierr = VecDestroy(&osm->ly);CHKERRQ(ierr); 698347574c9Seaulisa if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE) { 699347574c9Seaulisa //ierr = ISDestroy(&osm->lis);CHKERRQ(ierr); 700347574c9Seaulisa ierr = MatDestroyMatrices(osm->n_local_true, &osm->lmats);CHKERRQ(ierr); 701347574c9Seaulisa // ierr = VecDestroy(&osm->lx);CHKERRQ(ierr); 702347574c9Seaulisa // ierr = VecDestroy(&osm->ly);CHKERRQ(ierr); 703fb745f2cSMatthew G. Knepley } 7042fa5cd67SKarl Rupp 70580ec0b7dSPatrick Sanan ierr = PetscFree(osm->sub_mat_type);CHKERRQ(ierr); 70680ec0b7dSPatrick Sanan 707e91c6855SBarry Smith osm->is = 0; 708e91c6855SBarry Smith osm->is_local = 0; 709e91c6855SBarry Smith PetscFunctionReturn(0); 710e91c6855SBarry Smith } 711e91c6855SBarry Smith 712e91c6855SBarry Smith static PetscErrorCode PCDestroy_ASM(PC pc) 713e91c6855SBarry Smith { 714e91c6855SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 715e91c6855SBarry Smith PetscErrorCode ierr; 716e91c6855SBarry Smith PetscInt i; 717e91c6855SBarry Smith 718e91c6855SBarry Smith PetscFunctionBegin; 719e91c6855SBarry Smith ierr = PCReset_ASM(pc);CHKERRQ(ierr); 720e91c6855SBarry Smith if (osm->ksp) { 721e91c6855SBarry Smith for (i=0; i<osm->n_local_true; i++) { 722fcfd50ebSBarry Smith ierr = KSPDestroy(&osm->ksp[i]);CHKERRQ(ierr); 723e91c6855SBarry Smith } 724e91c6855SBarry Smith ierr = PetscFree(osm->ksp);CHKERRQ(ierr); 725e91c6855SBarry Smith } 726e91c6855SBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 7274b9ad928SBarry Smith PetscFunctionReturn(0); 7284b9ad928SBarry Smith } 7294b9ad928SBarry Smith 7304416b707SBarry Smith static PetscErrorCode PCSetFromOptions_ASM(PetscOptionItems *PetscOptionsObject,PC pc) 7314b9ad928SBarry Smith { 7324b9ad928SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 7336849ba73SBarry Smith PetscErrorCode ierr; 7349dcbbd2bSBarry Smith PetscInt blocks,ovl; 735ace3abfcSBarry Smith PetscBool symset,flg; 73692bb6962SLisandro Dalcin PCASMType asmtype; 73712cd4985SMatthew G. Knepley PCCompositeType loctype; 73880ec0b7dSPatrick Sanan char sub_mat_type[256]; 7394b9ad928SBarry Smith 7404b9ad928SBarry Smith PetscFunctionBegin; 741bf108f30SBarry Smith /* set the type to symmetric if matrix is symmetric */ 74292bb6962SLisandro Dalcin if (!osm->type_set && pc->pmat) { 74392bb6962SLisandro Dalcin ierr = MatIsSymmetricKnown(pc->pmat,&symset,&flg);CHKERRQ(ierr); 7442fa5cd67SKarl Rupp if (symset && flg) osm->type = PC_ASM_BASIC; 745bf108f30SBarry Smith } 746e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"Additive Schwarz options");CHKERRQ(ierr); 747d709ea83SDmitry Karpeev ierr = PetscOptionsBool("-pc_asm_dm_subdomains","Use DMCreateDomainDecomposition() to define subdomains","PCASMSetDMSubdomains",osm->dm_subdomains,&osm->dm_subdomains,&flg);CHKERRQ(ierr); 74890d69ab7SBarry Smith ierr = PetscOptionsInt("-pc_asm_blocks","Number of subdomains","PCASMSetTotalSubdomains",osm->n,&blocks,&flg);CHKERRQ(ierr); 74965db9045SDmitry Karpeev if (flg) { 750f77a5242SKarl Rupp ierr = PCASMSetTotalSubdomains(pc,blocks,NULL,NULL);CHKERRQ(ierr); 751d709ea83SDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 75265db9045SDmitry Karpeev } 75390d69ab7SBarry Smith ierr = PetscOptionsInt("-pc_asm_overlap","Number of grid points overlap","PCASMSetOverlap",osm->overlap,&ovl,&flg);CHKERRQ(ierr); 75465db9045SDmitry Karpeev if (flg) { 75565db9045SDmitry Karpeev ierr = PCASMSetOverlap(pc,ovl);CHKERRQ(ierr); 756d709ea83SDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 75765db9045SDmitry Karpeev } 75890d69ab7SBarry Smith flg = PETSC_FALSE; 75990d69ab7SBarry Smith ierr = PetscOptionsEnum("-pc_asm_type","Type of restriction/extension","PCASMSetType",PCASMTypes,(PetscEnum)osm->type,(PetscEnum*)&asmtype,&flg);CHKERRQ(ierr); 76092bb6962SLisandro Dalcin if (flg) {ierr = PCASMSetType(pc,asmtype);CHKERRQ(ierr); } 76112cd4985SMatthew G. Knepley flg = PETSC_FALSE; 76212cd4985SMatthew G. Knepley ierr = PetscOptionsEnum("-pc_asm_local_type","Type of local solver composition","PCASMSetLocalType",PCCompositeTypes,(PetscEnum)osm->loctype,(PetscEnum*)&loctype,&flg);CHKERRQ(ierr); 76312cd4985SMatthew G. Knepley if (flg) {ierr = PCASMSetLocalType(pc,loctype);CHKERRQ(ierr); } 76480ec0b7dSPatrick Sanan ierr = PetscOptionsFList("-pc_asm_sub_mat_type","Subsolve Matrix Type","PCASMSetSubMatType",MatList,NULL,sub_mat_type,256,&flg);CHKERRQ(ierr); 76580ec0b7dSPatrick Sanan if(flg){ 766459726d8SSatish Balay ierr = PCASMSetSubMatType(pc,sub_mat_type);CHKERRQ(ierr); 76780ec0b7dSPatrick Sanan } 7684b9ad928SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 7694b9ad928SBarry Smith PetscFunctionReturn(0); 7704b9ad928SBarry Smith } 7714b9ad928SBarry Smith 7724b9ad928SBarry Smith /*------------------------------------------------------------------------------------*/ 7734b9ad928SBarry Smith 7741e6b0712SBarry Smith static PetscErrorCode PCASMSetLocalSubdomains_ASM(PC pc,PetscInt n,IS is[],IS is_local[]) 7754b9ad928SBarry Smith { 7764b9ad928SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 77792bb6962SLisandro Dalcin PetscErrorCode ierr; 77892bb6962SLisandro Dalcin PetscInt i; 7794b9ad928SBarry Smith 7804b9ad928SBarry Smith PetscFunctionBegin; 781e32f2f54SBarry Smith if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Each process must have 1 or more blocks, n = %D",n); 782ce94432eSBarry Smith if (pc->setupcalled && (n != osm->n_local_true || is)) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCASMSetLocalSubdomains() should be called before calling PCSetUp()."); 783e7e72b3dSBarry Smith 7844b9ad928SBarry Smith if (!pc->setupcalled) { 78592bb6962SLisandro Dalcin if (is) { 78692bb6962SLisandro Dalcin for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)is[i]);CHKERRQ(ierr);} 78792bb6962SLisandro Dalcin } 788832fc9a5SMatthew Knepley if (is_local) { 789832fc9a5SMatthew Knepley for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)is_local[i]);CHKERRQ(ierr);} 790832fc9a5SMatthew Knepley } 7912b691e39SMatthew Knepley ierr = PCASMDestroySubdomains(osm->n_local_true,osm->is,osm->is_local);CHKERRQ(ierr); 7922fa5cd67SKarl Rupp 7934b9ad928SBarry Smith osm->n_local_true = n; 79492bb6962SLisandro Dalcin osm->is = 0; 7952b691e39SMatthew Knepley osm->is_local = 0; 79692bb6962SLisandro Dalcin if (is) { 797785e854fSJed Brown ierr = PetscMalloc1(n,&osm->is);CHKERRQ(ierr); 7982fa5cd67SKarl Rupp for (i=0; i<n; i++) osm->is[i] = is[i]; 7993d03bb48SJed Brown /* Flag indicating that the user has set overlapping subdomains so PCASM should not increase their size. */ 8003d03bb48SJed Brown osm->overlap = -1; 80192bb6962SLisandro Dalcin } 8022b691e39SMatthew Knepley if (is_local) { 803785e854fSJed Brown ierr = PetscMalloc1(n,&osm->is_local);CHKERRQ(ierr); 8042fa5cd67SKarl Rupp for (i=0; i<n; i++) osm->is_local[i] = is_local[i]; 805a35b7b57SDmitry Karpeev if (!is) { 806785e854fSJed Brown ierr = PetscMalloc1(osm->n_local_true,&osm->is);CHKERRQ(ierr); 807a35b7b57SDmitry Karpeev for (i=0; i<osm->n_local_true; i++) { 808a35b7b57SDmitry Karpeev if (osm->overlap > 0) { /* With positive overlap, osm->is[i] will be modified */ 809a35b7b57SDmitry Karpeev ierr = ISDuplicate(osm->is_local[i],&osm->is[i]);CHKERRQ(ierr); 810a35b7b57SDmitry Karpeev ierr = ISCopy(osm->is_local[i],osm->is[i]);CHKERRQ(ierr); 811a35b7b57SDmitry Karpeev } else { 812a35b7b57SDmitry Karpeev ierr = PetscObjectReference((PetscObject)osm->is_local[i]);CHKERRQ(ierr); 813a35b7b57SDmitry Karpeev osm->is[i] = osm->is_local[i]; 814a35b7b57SDmitry Karpeev } 815a35b7b57SDmitry Karpeev } 816a35b7b57SDmitry Karpeev } 8172b691e39SMatthew Knepley } 8184b9ad928SBarry Smith } 8194b9ad928SBarry Smith PetscFunctionReturn(0); 8204b9ad928SBarry Smith } 8214b9ad928SBarry Smith 8221e6b0712SBarry Smith static PetscErrorCode PCASMSetTotalSubdomains_ASM(PC pc,PetscInt N,IS *is,IS *is_local) 8234b9ad928SBarry Smith { 8244b9ad928SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 8256849ba73SBarry Smith PetscErrorCode ierr; 82613f74950SBarry Smith PetscMPIInt rank,size; 82778904715SLisandro Dalcin PetscInt n; 8284b9ad928SBarry Smith 8294b9ad928SBarry Smith PetscFunctionBegin; 830ce94432eSBarry Smith if (N < 1) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Number of total blocks must be > 0, N = %D",N); 831ce94432eSBarry Smith if (is || is_local) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Use PCASMSetLocalSubdomains() to set specific index sets\n\they cannot be set globally yet."); 8324b9ad928SBarry Smith 8334b9ad928SBarry Smith /* 834880770e9SJed Brown Split the subdomains equally among all processors 8354b9ad928SBarry Smith */ 836ce94432eSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);CHKERRQ(ierr); 837ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr); 8384b9ad928SBarry Smith n = N/size + ((N % size) > rank); 839e32f2f54SBarry Smith if (!n) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Process %d must have at least one block: total processors %d total blocks %D",(int)rank,(int)size,N); 840e32f2f54SBarry Smith if (pc->setupcalled && n != osm->n_local_true) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"PCASMSetTotalSubdomains() should be called before PCSetUp()."); 8414b9ad928SBarry Smith if (!pc->setupcalled) { 8422b691e39SMatthew Knepley ierr = PCASMDestroySubdomains(osm->n_local_true,osm->is,osm->is_local);CHKERRQ(ierr); 8432fa5cd67SKarl Rupp 8444b9ad928SBarry Smith osm->n_local_true = n; 8454b9ad928SBarry Smith osm->is = 0; 8462b691e39SMatthew Knepley osm->is_local = 0; 8474b9ad928SBarry Smith } 8484b9ad928SBarry Smith PetscFunctionReturn(0); 8494b9ad928SBarry Smith } 8504b9ad928SBarry Smith 8511e6b0712SBarry Smith static PetscErrorCode PCASMSetOverlap_ASM(PC pc,PetscInt ovl) 8524b9ad928SBarry Smith { 85392bb6962SLisandro Dalcin PC_ASM *osm = (PC_ASM*)pc->data; 8544b9ad928SBarry Smith 8554b9ad928SBarry Smith PetscFunctionBegin; 856ce94432eSBarry Smith if (ovl < 0) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap value requested"); 857ce94432eSBarry Smith if (pc->setupcalled && ovl != osm->overlap) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCASMSetOverlap() should be called before PCSetUp()."); 8582fa5cd67SKarl Rupp if (!pc->setupcalled) osm->overlap = ovl; 8594b9ad928SBarry Smith PetscFunctionReturn(0); 8604b9ad928SBarry Smith } 8614b9ad928SBarry Smith 8621e6b0712SBarry Smith static PetscErrorCode PCASMSetType_ASM(PC pc,PCASMType type) 8634b9ad928SBarry Smith { 86492bb6962SLisandro Dalcin PC_ASM *osm = (PC_ASM*)pc->data; 8654b9ad928SBarry Smith 8664b9ad928SBarry Smith PetscFunctionBegin; 8674b9ad928SBarry Smith osm->type = type; 868bf108f30SBarry Smith osm->type_set = PETSC_TRUE; 8694b9ad928SBarry Smith PetscFunctionReturn(0); 8704b9ad928SBarry Smith } 8714b9ad928SBarry Smith 872c60c7ad4SBarry Smith static PetscErrorCode PCASMGetType_ASM(PC pc,PCASMType *type) 873c60c7ad4SBarry Smith { 874c60c7ad4SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 875c60c7ad4SBarry Smith 876c60c7ad4SBarry Smith PetscFunctionBegin; 877c60c7ad4SBarry Smith *type = osm->type; 878c60c7ad4SBarry Smith PetscFunctionReturn(0); 879c60c7ad4SBarry Smith } 880c60c7ad4SBarry Smith 88112cd4985SMatthew G. Knepley static PetscErrorCode PCASMSetLocalType_ASM(PC pc, PCCompositeType type) 88212cd4985SMatthew G. Knepley { 88312cd4985SMatthew G. Knepley PC_ASM *osm = (PC_ASM *) pc->data; 88412cd4985SMatthew G. Knepley 88512cd4985SMatthew G. Knepley PetscFunctionBegin; 886d0ecd4dfSBarry Smith if (type != PC_COMPOSITE_ADDITIVE && type != PC_COMPOSITE_MULTIPLICATIVE) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Only supports additive or multiplicative as the local type"); 88712cd4985SMatthew G. Knepley osm->loctype = type; 88812cd4985SMatthew G. Knepley PetscFunctionReturn(0); 88912cd4985SMatthew G. Knepley } 89012cd4985SMatthew G. Knepley 89112cd4985SMatthew G. Knepley static PetscErrorCode PCASMGetLocalType_ASM(PC pc, PCCompositeType *type) 89212cd4985SMatthew G. Knepley { 89312cd4985SMatthew G. Knepley PC_ASM *osm = (PC_ASM *) pc->data; 89412cd4985SMatthew G. Knepley 89512cd4985SMatthew G. Knepley PetscFunctionBegin; 89612cd4985SMatthew G. Knepley *type = osm->loctype; 89712cd4985SMatthew G. Knepley PetscFunctionReturn(0); 89812cd4985SMatthew G. Knepley } 89912cd4985SMatthew G. Knepley 9001e6b0712SBarry Smith static PetscErrorCode PCASMSetSortIndices_ASM(PC pc,PetscBool doSort) 9016ed231c7SMatthew Knepley { 9026ed231c7SMatthew Knepley PC_ASM *osm = (PC_ASM*)pc->data; 9036ed231c7SMatthew Knepley 9046ed231c7SMatthew Knepley PetscFunctionBegin; 9056ed231c7SMatthew Knepley osm->sort_indices = doSort; 9066ed231c7SMatthew Knepley PetscFunctionReturn(0); 9076ed231c7SMatthew Knepley } 9086ed231c7SMatthew Knepley 9091e6b0712SBarry Smith static PetscErrorCode PCASMGetSubKSP_ASM(PC pc,PetscInt *n_local,PetscInt *first_local,KSP **ksp) 9104b9ad928SBarry Smith { 91192bb6962SLisandro Dalcin PC_ASM *osm = (PC_ASM*)pc->data; 912dfbe8321SBarry Smith PetscErrorCode ierr; 9134b9ad928SBarry Smith 9144b9ad928SBarry Smith PetscFunctionBegin; 915ce94432eSBarry Smith if (osm->n_local_true < 1) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Need to call PCSetUP() on PC (or KSPSetUp() on the outer KSP object) before calling here"); 9164b9ad928SBarry Smith 9172fa5cd67SKarl Rupp if (n_local) *n_local = osm->n_local_true; 91892bb6962SLisandro Dalcin if (first_local) { 919ce94432eSBarry Smith ierr = MPI_Scan(&osm->n_local_true,first_local,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr); 92092bb6962SLisandro Dalcin *first_local -= osm->n_local_true; 92192bb6962SLisandro Dalcin } 92292bb6962SLisandro Dalcin if (ksp) { 92392bb6962SLisandro Dalcin /* Assume that local solves are now different; not necessarily 92492bb6962SLisandro Dalcin true though! This flag is used only for PCView_ASM() */ 92592bb6962SLisandro Dalcin *ksp = osm->ksp; 92692bb6962SLisandro Dalcin osm->same_local_solves = PETSC_FALSE; 92792bb6962SLisandro Dalcin } 9284b9ad928SBarry Smith PetscFunctionReturn(0); 9294b9ad928SBarry Smith } 9304b9ad928SBarry Smith 93180ec0b7dSPatrick Sanan static PetscErrorCode PCASMGetSubMatType_ASM(PC pc,MatType *sub_mat_type) 93280ec0b7dSPatrick Sanan { 93380ec0b7dSPatrick Sanan PC_ASM *osm = (PC_ASM*)pc->data; 93480ec0b7dSPatrick Sanan 93580ec0b7dSPatrick Sanan PetscFunctionBegin; 93680ec0b7dSPatrick Sanan PetscValidHeaderSpecific(pc,PC_CLASSID,1); 93780ec0b7dSPatrick Sanan PetscValidPointer(sub_mat_type,2); 93880ec0b7dSPatrick Sanan *sub_mat_type = osm->sub_mat_type; 93980ec0b7dSPatrick Sanan PetscFunctionReturn(0); 94080ec0b7dSPatrick Sanan } 94180ec0b7dSPatrick Sanan 94280ec0b7dSPatrick Sanan static PetscErrorCode PCASMSetSubMatType_ASM(PC pc,MatType sub_mat_type) 94380ec0b7dSPatrick Sanan { 94480ec0b7dSPatrick Sanan PetscErrorCode ierr; 94580ec0b7dSPatrick Sanan PC_ASM *osm = (PC_ASM*)pc->data; 94680ec0b7dSPatrick Sanan 94780ec0b7dSPatrick Sanan PetscFunctionBegin; 94880ec0b7dSPatrick Sanan PetscValidHeaderSpecific(pc,PC_CLASSID,1); 94980ec0b7dSPatrick Sanan ierr = PetscFree(osm->sub_mat_type);CHKERRQ(ierr); 95080ec0b7dSPatrick Sanan ierr = PetscStrallocpy(sub_mat_type,(char**)&osm->sub_mat_type);CHKERRQ(ierr); 95180ec0b7dSPatrick Sanan PetscFunctionReturn(0); 95280ec0b7dSPatrick Sanan } 95380ec0b7dSPatrick Sanan 9544b9ad928SBarry Smith /*@C 9551093a601SBarry Smith PCASMSetLocalSubdomains - Sets the local subdomains (for this processor only) for the additive Schwarz preconditioner. 9564b9ad928SBarry Smith 9574b9ad928SBarry Smith Collective on PC 9584b9ad928SBarry Smith 9594b9ad928SBarry Smith Input Parameters: 9604b9ad928SBarry Smith + pc - the preconditioner context 9614b9ad928SBarry Smith . n - the number of subdomains for this processor (default value = 1) 9628c03b21aSDmitry Karpeev . is - the index set that defines the subdomains for this processor 9630298fd71SBarry Smith (or NULL for PETSc to determine subdomains) 964f1ee410cSBarry 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 965f1ee410cSBarry Smith (or NULL to not provide these) 9664b9ad928SBarry Smith 9674b9ad928SBarry Smith Notes: 9681093a601SBarry Smith The IS numbering is in the parallel, global numbering of the vector for both is and is_local 9694b9ad928SBarry Smith 9704b9ad928SBarry Smith By default the ASM preconditioner uses 1 block per processor. 9714b9ad928SBarry Smith 9724b9ad928SBarry Smith Use PCASMSetTotalSubdomains() to set the subdomains for all processors. 9734b9ad928SBarry Smith 974f1ee410cSBarry Smith If is_local is provided and PCASMType is PC_ASM_RESTRICT then the solution only over the is_local region is interpolated 975f1ee410cSBarry Smith back to form the global solution (this is the standard restricted additive Schwarz method) 976f1ee410cSBarry Smith 977f1ee410cSBarry 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 978f1ee410cSBarry Smith no code to handle that case. 979f1ee410cSBarry Smith 9804b9ad928SBarry Smith Level: advanced 9814b9ad928SBarry Smith 9824b9ad928SBarry Smith .keywords: PC, ASM, set, local, subdomains, additive Schwarz 9834b9ad928SBarry Smith 9844b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(), 985f1ee410cSBarry Smith PCASMCreateSubdomains2D(), PCASMGetLocalSubdomains(), PCASMType, PCASMSetType() 9864b9ad928SBarry Smith @*/ 9877087cfbeSBarry Smith PetscErrorCode PCASMSetLocalSubdomains(PC pc,PetscInt n,IS is[],IS is_local[]) 9884b9ad928SBarry Smith { 9897bb14e67SBarry Smith PetscErrorCode ierr; 9904b9ad928SBarry Smith 9914b9ad928SBarry Smith PetscFunctionBegin; 9920700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 9937bb14e67SBarry Smith ierr = PetscTryMethod(pc,"PCASMSetLocalSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,n,is,is_local));CHKERRQ(ierr); 9944b9ad928SBarry Smith PetscFunctionReturn(0); 9954b9ad928SBarry Smith } 9964b9ad928SBarry Smith 9974b9ad928SBarry Smith /*@C 998feb221f8SDmitry Karpeev PCASMSetTotalSubdomains - Sets the subdomains for all processors for the 9994b9ad928SBarry Smith additive Schwarz preconditioner. Either all or no processors in the 10004b9ad928SBarry Smith PC communicator must call this routine, with the same index sets. 10014b9ad928SBarry Smith 10024b9ad928SBarry Smith Collective on PC 10034b9ad928SBarry Smith 10044b9ad928SBarry Smith Input Parameters: 10054b9ad928SBarry Smith + pc - the preconditioner context 1006feb221f8SDmitry Karpeev . N - the number of subdomains for all processors 1007feb221f8SDmitry Karpeev . is - the index sets that define the subdomains for all processors 1008dfaa0c5fSBarry Smith (or NULL to ask PETSc to determine the subdomains) 10092b691e39SMatthew Knepley - is_local - the index sets that define the local part of the subdomains for this processor 1010f1ee410cSBarry Smith (or NULL to not provide this information) 10114b9ad928SBarry Smith 10124b9ad928SBarry Smith Options Database Key: 10134b9ad928SBarry Smith To set the total number of subdomain blocks rather than specify the 10144b9ad928SBarry Smith index sets, use the option 10154b9ad928SBarry Smith . -pc_asm_blocks <blks> - Sets total blocks 10164b9ad928SBarry Smith 10174b9ad928SBarry Smith Notes: 1018f1ee410cSBarry Smith Currently you cannot use this to set the actual subdomains with the argument is or is_local. 10194b9ad928SBarry Smith 10204b9ad928SBarry Smith By default the ASM preconditioner uses 1 block per processor. 10214b9ad928SBarry Smith 10224b9ad928SBarry Smith These index sets cannot be destroyed until after completion of the 10234b9ad928SBarry Smith linear solves for which the ASM preconditioner is being used. 10244b9ad928SBarry Smith 10254b9ad928SBarry Smith Use PCASMSetLocalSubdomains() to set local subdomains. 10264b9ad928SBarry Smith 10271093a601SBarry Smith The IS numbering is in the parallel, global numbering of the vector for both is and is_local 10281093a601SBarry Smith 10294b9ad928SBarry Smith Level: advanced 10304b9ad928SBarry Smith 10314b9ad928SBarry Smith .keywords: PC, ASM, set, total, global, subdomains, additive Schwarz 10324b9ad928SBarry Smith 10334b9ad928SBarry Smith .seealso: PCASMSetLocalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(), 10344b9ad928SBarry Smith PCASMCreateSubdomains2D() 10354b9ad928SBarry Smith @*/ 10367087cfbeSBarry Smith PetscErrorCode PCASMSetTotalSubdomains(PC pc,PetscInt N,IS is[],IS is_local[]) 10374b9ad928SBarry Smith { 10387bb14e67SBarry Smith PetscErrorCode ierr; 10394b9ad928SBarry Smith 10404b9ad928SBarry Smith PetscFunctionBegin; 10410700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 10427bb14e67SBarry Smith ierr = PetscTryMethod(pc,"PCASMSetTotalSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,N,is,is_local));CHKERRQ(ierr); 10434b9ad928SBarry Smith PetscFunctionReturn(0); 10444b9ad928SBarry Smith } 10454b9ad928SBarry Smith 10464b9ad928SBarry Smith /*@ 10474b9ad928SBarry Smith PCASMSetOverlap - Sets the overlap between a pair of subdomains for the 10484b9ad928SBarry Smith additive Schwarz preconditioner. Either all or no processors in the 1049f1ee410cSBarry Smith PC communicator must call this routine. 10504b9ad928SBarry Smith 1051ad4df100SBarry Smith Logically Collective on PC 10524b9ad928SBarry Smith 10534b9ad928SBarry Smith Input Parameters: 10544b9ad928SBarry Smith + pc - the preconditioner context 10554b9ad928SBarry Smith - ovl - the amount of overlap between subdomains (ovl >= 0, default value = 1) 10564b9ad928SBarry Smith 10574b9ad928SBarry Smith Options Database Key: 10584b9ad928SBarry Smith . -pc_asm_overlap <ovl> - Sets overlap 10594b9ad928SBarry Smith 10604b9ad928SBarry Smith Notes: 10614b9ad928SBarry Smith By default the ASM preconditioner uses 1 block per processor. To use 10624b9ad928SBarry Smith multiple blocks per perocessor, see PCASMSetTotalSubdomains() and 10634b9ad928SBarry Smith PCASMSetLocalSubdomains() (and the option -pc_asm_blocks <blks>). 10644b9ad928SBarry Smith 10654b9ad928SBarry Smith The overlap defaults to 1, so if one desires that no additional 10664b9ad928SBarry Smith overlap be computed beyond what may have been set with a call to 10674b9ad928SBarry Smith PCASMSetTotalSubdomains() or PCASMSetLocalSubdomains(), then ovl 10684b9ad928SBarry Smith must be set to be 0. In particular, if one does not explicitly set 10694b9ad928SBarry Smith the subdomains an application code, then all overlap would be computed 10704b9ad928SBarry Smith internally by PETSc, and using an overlap of 0 would result in an ASM 10714b9ad928SBarry Smith variant that is equivalent to the block Jacobi preconditioner. 10724b9ad928SBarry Smith 1073f1ee410cSBarry Smith The default algorithm used by PETSc to increase overlap is fast, but not scalable, 1074f1ee410cSBarry Smith use the option -mat_increase_overlap_scalable when the problem and number of processes is large. 1075f1ee410cSBarry Smith 10764b9ad928SBarry Smith Note that one can define initial index sets with any overlap via 1077f1ee410cSBarry Smith PCASMSetLocalSubdomains(); the routine 10784b9ad928SBarry Smith PCASMSetOverlap() merely allows PETSc to extend that overlap further 10794b9ad928SBarry Smith if desired. 10804b9ad928SBarry Smith 10814b9ad928SBarry Smith Level: intermediate 10824b9ad928SBarry Smith 10834b9ad928SBarry Smith .keywords: PC, ASM, set, overlap 10844b9ad928SBarry Smith 10854b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetLocalSubdomains(), PCASMGetSubKSP(), 1086f1ee410cSBarry Smith PCASMCreateSubdomains2D(), PCASMGetLocalSubdomains(), MatIncreaseOverlap() 10874b9ad928SBarry Smith @*/ 10887087cfbeSBarry Smith PetscErrorCode PCASMSetOverlap(PC pc,PetscInt ovl) 10894b9ad928SBarry Smith { 10907bb14e67SBarry Smith PetscErrorCode ierr; 10914b9ad928SBarry Smith 10924b9ad928SBarry Smith PetscFunctionBegin; 10930700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1094c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,ovl,2); 10957bb14e67SBarry Smith ierr = PetscTryMethod(pc,"PCASMSetOverlap_C",(PC,PetscInt),(pc,ovl));CHKERRQ(ierr); 10964b9ad928SBarry Smith PetscFunctionReturn(0); 10974b9ad928SBarry Smith } 10984b9ad928SBarry Smith 10994b9ad928SBarry Smith /*@ 11004b9ad928SBarry Smith PCASMSetType - Sets the type of restriction and interpolation used 11014b9ad928SBarry Smith for local problems in the additive Schwarz method. 11024b9ad928SBarry Smith 1103ad4df100SBarry Smith Logically Collective on PC 11044b9ad928SBarry Smith 11054b9ad928SBarry Smith Input Parameters: 11064b9ad928SBarry Smith + pc - the preconditioner context 11074b9ad928SBarry Smith - type - variant of ASM, one of 11084b9ad928SBarry Smith .vb 11094b9ad928SBarry Smith PC_ASM_BASIC - full interpolation and restriction 11104b9ad928SBarry Smith PC_ASM_RESTRICT - full restriction, local processor interpolation 11114b9ad928SBarry Smith PC_ASM_INTERPOLATE - full interpolation, local processor restriction 11124b9ad928SBarry Smith PC_ASM_NONE - local processor restriction and interpolation 11134b9ad928SBarry Smith .ve 11144b9ad928SBarry Smith 11154b9ad928SBarry Smith Options Database Key: 11164b9ad928SBarry Smith . -pc_asm_type [basic,restrict,interpolate,none] - Sets ASM type 11174b9ad928SBarry Smith 1118f1ee410cSBarry Smith Notes: if the is_local arguments are passed to PCASMSetLocalSubdomains() then they are used when PC_ASM_RESTRICT has been selected 1119f1ee410cSBarry Smith to limit the local processor interpolation 1120f1ee410cSBarry Smith 11214b9ad928SBarry Smith Level: intermediate 11224b9ad928SBarry Smith 11234b9ad928SBarry Smith .keywords: PC, ASM, set, type 11244b9ad928SBarry Smith 11254b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetTotalSubdomains(), PCASMGetSubKSP(), 1126f1ee410cSBarry Smith PCASMCreateSubdomains2D(), PCASMType, PCASMSetLocalType(), PCASMGetLocalType() 11274b9ad928SBarry Smith @*/ 11287087cfbeSBarry Smith PetscErrorCode PCASMSetType(PC pc,PCASMType type) 11294b9ad928SBarry Smith { 11307bb14e67SBarry Smith PetscErrorCode ierr; 11314b9ad928SBarry Smith 11324b9ad928SBarry Smith PetscFunctionBegin; 11330700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1134c5eb9154SBarry Smith PetscValidLogicalCollectiveEnum(pc,type,2); 11357bb14e67SBarry Smith ierr = PetscTryMethod(pc,"PCASMSetType_C",(PC,PCASMType),(pc,type));CHKERRQ(ierr); 11364b9ad928SBarry Smith PetscFunctionReturn(0); 11374b9ad928SBarry Smith } 11384b9ad928SBarry Smith 1139c60c7ad4SBarry Smith /*@ 1140c60c7ad4SBarry Smith PCASMGetType - Gets the type of restriction and interpolation used 1141c60c7ad4SBarry Smith for local problems in the additive Schwarz method. 1142c60c7ad4SBarry Smith 1143c60c7ad4SBarry Smith Logically Collective on PC 1144c60c7ad4SBarry Smith 1145c60c7ad4SBarry Smith Input Parameter: 1146c60c7ad4SBarry Smith . pc - the preconditioner context 1147c60c7ad4SBarry Smith 1148c60c7ad4SBarry Smith Output Parameter: 1149c60c7ad4SBarry Smith . type - variant of ASM, one of 1150c60c7ad4SBarry Smith 1151c60c7ad4SBarry Smith .vb 1152c60c7ad4SBarry Smith PC_ASM_BASIC - full interpolation and restriction 1153c60c7ad4SBarry Smith PC_ASM_RESTRICT - full restriction, local processor interpolation 1154c60c7ad4SBarry Smith PC_ASM_INTERPOLATE - full interpolation, local processor restriction 1155c60c7ad4SBarry Smith PC_ASM_NONE - local processor restriction and interpolation 1156c60c7ad4SBarry Smith .ve 1157c60c7ad4SBarry Smith 1158c60c7ad4SBarry Smith Options Database Key: 1159c60c7ad4SBarry Smith . -pc_asm_type [basic,restrict,interpolate,none] - Sets ASM type 1160c60c7ad4SBarry Smith 1161c60c7ad4SBarry Smith Level: intermediate 1162c60c7ad4SBarry Smith 1163c60c7ad4SBarry Smith .keywords: PC, ASM, set, type 1164c60c7ad4SBarry Smith 1165c60c7ad4SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetTotalSubdomains(), PCASMGetSubKSP(), 1166f1ee410cSBarry Smith PCASMCreateSubdomains2D(), PCASMType, PCASMSetType(), PCASMSetLocalType(), PCASMGetLocalType() 1167c60c7ad4SBarry Smith @*/ 1168c60c7ad4SBarry Smith PetscErrorCode PCASMGetType(PC pc,PCASMType *type) 1169c60c7ad4SBarry Smith { 1170c60c7ad4SBarry Smith PetscErrorCode ierr; 1171c60c7ad4SBarry Smith 1172c60c7ad4SBarry Smith PetscFunctionBegin; 1173c60c7ad4SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1174c60c7ad4SBarry Smith ierr = PetscUseMethod(pc,"PCASMGetType_C",(PC,PCASMType*),(pc,type));CHKERRQ(ierr); 1175c60c7ad4SBarry Smith PetscFunctionReturn(0); 1176c60c7ad4SBarry Smith } 1177c60c7ad4SBarry Smith 117812cd4985SMatthew G. Knepley /*@ 117912cd4985SMatthew G. Knepley PCASMSetLocalType - Sets the type of composition used for local problems in the additive Schwarz method. 118012cd4985SMatthew G. Knepley 118112cd4985SMatthew G. Knepley Logically Collective on PC 118212cd4985SMatthew G. Knepley 118312cd4985SMatthew G. Knepley Input Parameters: 118412cd4985SMatthew G. Knepley + pc - the preconditioner context 118512cd4985SMatthew G. Knepley - type - type of composition, one of 118612cd4985SMatthew G. Knepley .vb 118712cd4985SMatthew G. Knepley PC_COMPOSITE_ADDITIVE - local additive combination 118812cd4985SMatthew G. Knepley PC_COMPOSITE_MULTIPLICATIVE - local multiplicative combination 118912cd4985SMatthew G. Knepley .ve 119012cd4985SMatthew G. Knepley 119112cd4985SMatthew G. Knepley Options Database Key: 119212cd4985SMatthew G. Knepley . -pc_asm_local_type [additive,multiplicative] - Sets local solver composition type 119312cd4985SMatthew G. Knepley 119412cd4985SMatthew G. Knepley Level: intermediate 119512cd4985SMatthew G. Knepley 1196f1ee410cSBarry Smith .seealso: PCASMSetType(), PCASMGetType(), PCASMGetLocalType(), PCASM, PCASMType, PCASMSetType(), PCASMGetType(), PCCompositeType 119712cd4985SMatthew G. Knepley @*/ 119812cd4985SMatthew G. Knepley PetscErrorCode PCASMSetLocalType(PC pc, PCCompositeType type) 119912cd4985SMatthew G. Knepley { 120012cd4985SMatthew G. Knepley PetscErrorCode ierr; 120112cd4985SMatthew G. Knepley 120212cd4985SMatthew G. Knepley PetscFunctionBegin; 120312cd4985SMatthew G. Knepley PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 120412cd4985SMatthew G. Knepley PetscValidLogicalCollectiveEnum(pc, type, 2); 120512cd4985SMatthew G. Knepley ierr = PetscTryMethod(pc, "PCASMSetLocalType_C", (PC, PCCompositeType), (pc, type));CHKERRQ(ierr); 120612cd4985SMatthew G. Knepley PetscFunctionReturn(0); 120712cd4985SMatthew G. Knepley } 120812cd4985SMatthew G. Knepley 120912cd4985SMatthew G. Knepley /*@ 121012cd4985SMatthew G. Knepley PCASMGetLocalType - Gets the type of composition used for local problems in the additive Schwarz method. 121112cd4985SMatthew G. Knepley 121212cd4985SMatthew G. Knepley Logically Collective on PC 121312cd4985SMatthew G. Knepley 121412cd4985SMatthew G. Knepley Input Parameter: 121512cd4985SMatthew G. Knepley . pc - the preconditioner context 121612cd4985SMatthew G. Knepley 121712cd4985SMatthew G. Knepley Output Parameter: 121812cd4985SMatthew G. Knepley . type - type of composition, one of 121912cd4985SMatthew G. Knepley .vb 122012cd4985SMatthew G. Knepley PC_COMPOSITE_ADDITIVE - local additive combination 122112cd4985SMatthew G. Knepley PC_COMPOSITE_MULTIPLICATIVE - local multiplicative combination 122212cd4985SMatthew G. Knepley .ve 122312cd4985SMatthew G. Knepley 122412cd4985SMatthew G. Knepley Options Database Key: 122512cd4985SMatthew G. Knepley . -pc_asm_local_type [additive,multiplicative] - Sets local solver composition type 122612cd4985SMatthew G. Knepley 122712cd4985SMatthew G. Knepley Level: intermediate 122812cd4985SMatthew G. Knepley 1229f1ee410cSBarry Smith .seealso: PCASMSetType(), PCASMGetType(), PCASMSetLocalType(), PCASMCreate(), PCASMType, PCASMSetType(), PCASMGetType(), PCCompositeType 123012cd4985SMatthew G. Knepley @*/ 123112cd4985SMatthew G. Knepley PetscErrorCode PCASMGetLocalType(PC pc, PCCompositeType *type) 123212cd4985SMatthew G. Knepley { 123312cd4985SMatthew G. Knepley PetscErrorCode ierr; 123412cd4985SMatthew G. Knepley 123512cd4985SMatthew G. Knepley PetscFunctionBegin; 123612cd4985SMatthew G. Knepley PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 123712cd4985SMatthew G. Knepley PetscValidPointer(type, 2); 123812cd4985SMatthew G. Knepley ierr = PetscUseMethod(pc, "PCASMGetLocalType_C", (PC, PCCompositeType *), (pc, type));CHKERRQ(ierr); 123912cd4985SMatthew G. Knepley PetscFunctionReturn(0); 124012cd4985SMatthew G. Knepley } 124112cd4985SMatthew G. Knepley 12426ed231c7SMatthew Knepley /*@ 12436ed231c7SMatthew Knepley PCASMSetSortIndices - Determines whether subdomain indices are sorted. 12446ed231c7SMatthew Knepley 1245ad4df100SBarry Smith Logically Collective on PC 12466ed231c7SMatthew Knepley 12476ed231c7SMatthew Knepley Input Parameters: 12486ed231c7SMatthew Knepley + pc - the preconditioner context 12496ed231c7SMatthew Knepley - doSort - sort the subdomain indices 12506ed231c7SMatthew Knepley 12516ed231c7SMatthew Knepley Level: intermediate 12526ed231c7SMatthew Knepley 12536ed231c7SMatthew Knepley .keywords: PC, ASM, set, type 12546ed231c7SMatthew Knepley 12556ed231c7SMatthew Knepley .seealso: PCASMSetLocalSubdomains(), PCASMSetTotalSubdomains(), PCASMGetSubKSP(), 12566ed231c7SMatthew Knepley PCASMCreateSubdomains2D() 12576ed231c7SMatthew Knepley @*/ 12587087cfbeSBarry Smith PetscErrorCode PCASMSetSortIndices(PC pc,PetscBool doSort) 12596ed231c7SMatthew Knepley { 12607bb14e67SBarry Smith PetscErrorCode ierr; 12616ed231c7SMatthew Knepley 12626ed231c7SMatthew Knepley PetscFunctionBegin; 12630700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1264acfcf0e5SJed Brown PetscValidLogicalCollectiveBool(pc,doSort,2); 12657bb14e67SBarry Smith ierr = PetscTryMethod(pc,"PCASMSetSortIndices_C",(PC,PetscBool),(pc,doSort));CHKERRQ(ierr); 12666ed231c7SMatthew Knepley PetscFunctionReturn(0); 12676ed231c7SMatthew Knepley } 12686ed231c7SMatthew Knepley 12694b9ad928SBarry Smith /*@C 12704b9ad928SBarry Smith PCASMGetSubKSP - Gets the local KSP contexts for all blocks on 12714b9ad928SBarry Smith this processor. 12724b9ad928SBarry Smith 12734b9ad928SBarry Smith Collective on PC iff first_local is requested 12744b9ad928SBarry Smith 12754b9ad928SBarry Smith Input Parameter: 12764b9ad928SBarry Smith . pc - the preconditioner context 12774b9ad928SBarry Smith 12784b9ad928SBarry Smith Output Parameters: 12790298fd71SBarry Smith + n_local - the number of blocks on this processor or NULL 12800298fd71SBarry Smith . first_local - the global number of the first block on this processor or NULL, 12810298fd71SBarry Smith all processors must request or all must pass NULL 12824b9ad928SBarry Smith - ksp - the array of KSP contexts 12834b9ad928SBarry Smith 12844b9ad928SBarry Smith Note: 1285d29017ddSJed Brown After PCASMGetSubKSP() the array of KSPes is not to be freed. 12864b9ad928SBarry Smith 12874b9ad928SBarry Smith You must call KSPSetUp() before calling PCASMGetSubKSP(). 12884b9ad928SBarry Smith 1289d29017ddSJed Brown Fortran note: 12902bf68e3eSBarry 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. 1291d29017ddSJed Brown 12924b9ad928SBarry Smith Level: advanced 12934b9ad928SBarry Smith 12944b9ad928SBarry Smith .keywords: PC, ASM, additive Schwarz, get, sub, KSP, context 12954b9ad928SBarry Smith 12964b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetTotalSubdomains(), PCASMSetOverlap(), 12974b9ad928SBarry Smith PCASMCreateSubdomains2D(), 12984b9ad928SBarry Smith @*/ 12997087cfbeSBarry Smith PetscErrorCode PCASMGetSubKSP(PC pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp[]) 13004b9ad928SBarry Smith { 13017bb14e67SBarry Smith PetscErrorCode ierr; 13024b9ad928SBarry Smith 13034b9ad928SBarry Smith PetscFunctionBegin; 13040700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 13057bb14e67SBarry Smith ierr = PetscUseMethod(pc,"PCASMGetSubKSP_C",(PC,PetscInt*,PetscInt*,KSP **),(pc,n_local,first_local,ksp));CHKERRQ(ierr); 13064b9ad928SBarry Smith PetscFunctionReturn(0); 13074b9ad928SBarry Smith } 13084b9ad928SBarry Smith 13094b9ad928SBarry Smith /* -------------------------------------------------------------------------------------*/ 13104b9ad928SBarry Smith /*MC 13113b09bd56SBarry Smith PCASM - Use the (restricted) additive Schwarz method, each block is (approximately) solved with 13124b9ad928SBarry Smith its own KSP object. 13134b9ad928SBarry Smith 13144b9ad928SBarry Smith Options Database Keys: 131549517cdeSBarry Smith + -pc_asm_blocks <blks> - Sets total blocks 13164b9ad928SBarry Smith . -pc_asm_overlap <ovl> - Sets overlap 1317f1ee410cSBarry Smith . -pc_asm_type [basic,restrict,interpolate,none] - Sets ASM type, default is restrict 1318f1ee410cSBarry Smith - -pc_asm_local_type [additive, multiplicative] - Sets ASM type, default is additive 13194b9ad928SBarry Smith 13203b09bd56SBarry Smith IMPORTANT: If you run with, for example, 3 blocks on 1 processor or 3 blocks on 3 processors you 13213b09bd56SBarry Smith will get a different convergence rate due to the default option of -pc_asm_type restrict. Use 13223b09bd56SBarry Smith -pc_asm_type basic to use the standard ASM. 13233b09bd56SBarry Smith 13244b9ad928SBarry Smith Notes: Each processor can have one or more blocks, but a block cannot be shared by more 1325f1ee410cSBarry Smith than one processor. Use PCGASM for subdomains shared by multiple processes. Defaults to one block per processor. 13264b9ad928SBarry Smith 13273b09bd56SBarry Smith To set options on the solvers for each block append -sub_ to all the KSP, and PC 1328d7ee0231SBarry Smith options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly 13294b9ad928SBarry Smith 1330a8c7a070SBarry Smith To set the options on the solvers separate for each block call PCASMGetSubKSP() 13314b9ad928SBarry Smith and set the options directly on the resulting KSP object (you can access its PC 13324b9ad928SBarry Smith with KSPGetPC()) 13334b9ad928SBarry Smith 13344b9ad928SBarry Smith Level: beginner 13354b9ad928SBarry Smith 13364b9ad928SBarry Smith Concepts: additive Schwarz method 13374b9ad928SBarry Smith 1338c582cd25SBarry Smith References: 133996a0c994SBarry Smith + 1. - M Dryja, OB Widlund, An additive variant of the Schwarz alternating method for the case of many subregions 134096a0c994SBarry Smith Courant Institute, New York University Technical report 134196a0c994SBarry Smith - 1. - Barry Smith, Petter Bjorstad, and William Gropp, Domain Decompositions: Parallel Multilevel Methods for Elliptic Partial Differential Equations, 134296a0c994SBarry Smith Cambridge University Press. 1343c582cd25SBarry Smith 13444b9ad928SBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 1345f1ee410cSBarry Smith PCBJACOBI, PCASMGetSubKSP(), PCASMSetLocalSubdomains(), PCASMType, PCASMGetType(), PCASMSetLocalType(), PCASMGetLocalType() 1346f1ee410cSBarry Smith PCASMSetTotalSubdomains(), PCSetModifySubmatrices(), PCASMSetOverlap(), PCASMSetType(), PCCompositeType 1347e09e08ccSBarry Smith 13484b9ad928SBarry Smith M*/ 13494b9ad928SBarry Smith 13508cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_ASM(PC pc) 13514b9ad928SBarry Smith { 1352dfbe8321SBarry Smith PetscErrorCode ierr; 13534b9ad928SBarry Smith PC_ASM *osm; 13544b9ad928SBarry Smith 13554b9ad928SBarry Smith PetscFunctionBegin; 1356b00a9115SJed Brown ierr = PetscNewLog(pc,&osm);CHKERRQ(ierr); 13572fa5cd67SKarl Rupp 13584b9ad928SBarry Smith osm->n = PETSC_DECIDE; 13594b9ad928SBarry Smith osm->n_local = 0; 13602b837212SDmitry Karpeev osm->n_local_true = PETSC_DECIDE; 13614b9ad928SBarry Smith osm->overlap = 1; 13624b9ad928SBarry Smith osm->ksp = 0; 13632b691e39SMatthew Knepley osm->restriction = 0; 1364347574c9Seaulisa osm->lrestriction = 0; 1365b9845e0eSMatthew Knepley osm->localization = 0; 13662b691e39SMatthew Knepley osm->prolongation = 0; 13675b3c0d42Seaulisa osm->lprolongation = 0; 136892bb6962SLisandro Dalcin osm->x = 0; 136992bb6962SLisandro Dalcin osm->y = 0; 137053888de8SMatthew Knepley osm->y_local = 0; 13714b9ad928SBarry Smith osm->is = 0; 13722b691e39SMatthew Knepley osm->is_local = 0; 13734b9ad928SBarry Smith osm->mat = 0; 13744b9ad928SBarry Smith osm->pmat = 0; 13754b9ad928SBarry Smith osm->type = PC_ASM_RESTRICT; 137612cd4985SMatthew G. Knepley osm->loctype = PC_COMPOSITE_ADDITIVE; 13774b9ad928SBarry Smith osm->same_local_solves = PETSC_TRUE; 13786ed231c7SMatthew Knepley osm->sort_indices = PETSC_TRUE; 1379d709ea83SDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 138080ec0b7dSPatrick Sanan osm->sub_mat_type = NULL; 13814b9ad928SBarry Smith 138292bb6962SLisandro Dalcin pc->data = (void*)osm; 13834b9ad928SBarry Smith pc->ops->apply = PCApply_ASM; 13844b9ad928SBarry Smith pc->ops->applytranspose = PCApplyTranspose_ASM; 13854b9ad928SBarry Smith pc->ops->setup = PCSetUp_ASM; 1386e91c6855SBarry Smith pc->ops->reset = PCReset_ASM; 13874b9ad928SBarry Smith pc->ops->destroy = PCDestroy_ASM; 13884b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_ASM; 13894b9ad928SBarry Smith pc->ops->setuponblocks = PCSetUpOnBlocks_ASM; 13904b9ad928SBarry Smith pc->ops->view = PCView_ASM; 13914b9ad928SBarry Smith pc->ops->applyrichardson = 0; 13924b9ad928SBarry Smith 1393bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetLocalSubdomains_C",PCASMSetLocalSubdomains_ASM);CHKERRQ(ierr); 1394bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetTotalSubdomains_C",PCASMSetTotalSubdomains_ASM);CHKERRQ(ierr); 1395bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetOverlap_C",PCASMSetOverlap_ASM);CHKERRQ(ierr); 1396bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetType_C",PCASMSetType_ASM);CHKERRQ(ierr); 1397c60c7ad4SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMGetType_C",PCASMGetType_ASM);CHKERRQ(ierr); 139812cd4985SMatthew G. Knepley ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetLocalType_C",PCASMSetLocalType_ASM);CHKERRQ(ierr); 139912cd4985SMatthew G. Knepley ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMGetLocalType_C",PCASMGetLocalType_ASM);CHKERRQ(ierr); 1400bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetSortIndices_C",PCASMSetSortIndices_ASM);CHKERRQ(ierr); 1401bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMGetSubKSP_C",PCASMGetSubKSP_ASM);CHKERRQ(ierr); 140280ec0b7dSPatrick Sanan ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMGetSubMatType_C",PCASMGetSubMatType_ASM);CHKERRQ(ierr); 140380ec0b7dSPatrick Sanan ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetSubMatType_C",PCASMSetSubMatType_ASM);CHKERRQ(ierr); 14044b9ad928SBarry Smith PetscFunctionReturn(0); 14054b9ad928SBarry Smith } 14064b9ad928SBarry Smith 140792bb6962SLisandro Dalcin /*@C 140892bb6962SLisandro Dalcin PCASMCreateSubdomains - Creates the index sets for the overlapping Schwarz 140992bb6962SLisandro Dalcin preconditioner for a any problem on a general grid. 141092bb6962SLisandro Dalcin 141192bb6962SLisandro Dalcin Collective 141292bb6962SLisandro Dalcin 141392bb6962SLisandro Dalcin Input Parameters: 141492bb6962SLisandro Dalcin + A - The global matrix operator 141592bb6962SLisandro Dalcin - n - the number of local blocks 141692bb6962SLisandro Dalcin 141792bb6962SLisandro Dalcin Output Parameters: 141892bb6962SLisandro Dalcin . outis - the array of index sets defining the subdomains 141992bb6962SLisandro Dalcin 142092bb6962SLisandro Dalcin Level: advanced 142192bb6962SLisandro Dalcin 14227d6bfa3bSBarry Smith Note: this generates nonoverlapping subdomains; the PCASM will generate the overlap 14237d6bfa3bSBarry Smith from these if you use PCASMSetLocalSubdomains() 14247d6bfa3bSBarry Smith 14257d6bfa3bSBarry Smith In the Fortran version you must provide the array outis[] already allocated of length n. 14267d6bfa3bSBarry Smith 142792bb6962SLisandro Dalcin .keywords: PC, ASM, additive Schwarz, create, subdomains, unstructured grid 142892bb6962SLisandro Dalcin 142992bb6962SLisandro Dalcin .seealso: PCASMSetLocalSubdomains(), PCASMDestroySubdomains() 143092bb6962SLisandro Dalcin @*/ 14317087cfbeSBarry Smith PetscErrorCode PCASMCreateSubdomains(Mat A, PetscInt n, IS* outis[]) 143292bb6962SLisandro Dalcin { 143392bb6962SLisandro Dalcin MatPartitioning mpart; 143492bb6962SLisandro Dalcin const char *prefix; 1435c56a70eeSBarry Smith void (*f)(void); 143692bb6962SLisandro Dalcin PetscInt i,j,rstart,rend,bs; 143711bd1e4dSLisandro Dalcin PetscBool isbaij = PETSC_FALSE,foundpart = PETSC_FALSE; 14380298fd71SBarry Smith Mat Ad = NULL, adj; 143992bb6962SLisandro Dalcin IS ispart,isnumb,*is; 144092bb6962SLisandro Dalcin PetscErrorCode ierr; 144192bb6962SLisandro Dalcin 144292bb6962SLisandro Dalcin PetscFunctionBegin; 14430700a824SBarry Smith PetscValidHeaderSpecific(A,MAT_CLASSID,1); 144492bb6962SLisandro Dalcin PetscValidPointer(outis,3); 1445c1235816SBarry Smith if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of local blocks must be > 0, n = %D",n); 144692bb6962SLisandro Dalcin 144792bb6962SLisandro Dalcin /* Get prefix, row distribution, and block size */ 144892bb6962SLisandro Dalcin ierr = MatGetOptionsPrefix(A,&prefix);CHKERRQ(ierr); 144992bb6962SLisandro Dalcin ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 145092bb6962SLisandro Dalcin ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr); 145165e19b50SBarry Smith if (rstart/bs*bs != rstart || rend/bs*bs != rend) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"bad row distribution [%D,%D) for matrix block size %D",rstart,rend,bs); 145265e19b50SBarry Smith 145392bb6962SLisandro Dalcin /* Get diagonal block from matrix if possible */ 1454c56a70eeSBarry Smith ierr = MatShellGetOperation(A,MATOP_GET_DIAGONAL_BLOCK,&f);CHKERRQ(ierr); 145592bb6962SLisandro Dalcin if (f) { 145611bd1e4dSLisandro Dalcin ierr = MatGetDiagonalBlock(A,&Ad);CHKERRQ(ierr); 145792bb6962SLisandro Dalcin } 145892bb6962SLisandro Dalcin if (Ad) { 1459251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)Ad,MATSEQBAIJ,&isbaij);CHKERRQ(ierr); 1460251f4c67SDmitry Karpeev if (!isbaij) {ierr = PetscObjectTypeCompare((PetscObject)Ad,MATSEQSBAIJ,&isbaij);CHKERRQ(ierr);} 146192bb6962SLisandro Dalcin } 146292bb6962SLisandro Dalcin if (Ad && n > 1) { 1463ace3abfcSBarry Smith PetscBool match,done; 146492bb6962SLisandro Dalcin /* Try to setup a good matrix partitioning if available */ 146592bb6962SLisandro Dalcin ierr = MatPartitioningCreate(PETSC_COMM_SELF,&mpart);CHKERRQ(ierr); 146692bb6962SLisandro Dalcin ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr); 146792bb6962SLisandro Dalcin ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr); 1468251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGCURRENT,&match);CHKERRQ(ierr); 146992bb6962SLisandro Dalcin if (!match) { 1470251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGSQUARE,&match);CHKERRQ(ierr); 147192bb6962SLisandro Dalcin } 147292bb6962SLisandro Dalcin if (!match) { /* assume a "good" partitioner is available */ 14731a83f524SJed Brown PetscInt na; 14741a83f524SJed Brown const PetscInt *ia,*ja; 147592bb6962SLisandro Dalcin ierr = MatGetRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr); 147692bb6962SLisandro Dalcin if (done) { 147792bb6962SLisandro Dalcin /* Build adjacency matrix by hand. Unfortunately a call to 147892bb6962SLisandro Dalcin MatConvert(Ad,MATMPIADJ,MAT_INITIAL_MATRIX,&adj) will 147992bb6962SLisandro Dalcin remove the block-aij structure and we cannot expect 148092bb6962SLisandro Dalcin MatPartitioning to split vertices as we need */ 14811a83f524SJed Brown PetscInt i,j,len,nnz,cnt,*iia=0,*jja=0; 14821a83f524SJed Brown const PetscInt *row; 148392bb6962SLisandro Dalcin nnz = 0; 148492bb6962SLisandro Dalcin for (i=0; i<na; i++) { /* count number of nonzeros */ 148592bb6962SLisandro Dalcin len = ia[i+1] - ia[i]; 148692bb6962SLisandro Dalcin row = ja + ia[i]; 148792bb6962SLisandro Dalcin for (j=0; j<len; j++) { 148892bb6962SLisandro Dalcin if (row[j] == i) { /* don't count diagonal */ 148992bb6962SLisandro Dalcin len--; break; 149092bb6962SLisandro Dalcin } 149192bb6962SLisandro Dalcin } 149292bb6962SLisandro Dalcin nnz += len; 149392bb6962SLisandro Dalcin } 1494854ce69bSBarry Smith ierr = PetscMalloc1(na+1,&iia);CHKERRQ(ierr); 1495854ce69bSBarry Smith ierr = PetscMalloc1(nnz,&jja);CHKERRQ(ierr); 149692bb6962SLisandro Dalcin nnz = 0; 149792bb6962SLisandro Dalcin iia[0] = 0; 149892bb6962SLisandro Dalcin for (i=0; i<na; i++) { /* fill adjacency */ 149992bb6962SLisandro Dalcin cnt = 0; 150092bb6962SLisandro Dalcin len = ia[i+1] - ia[i]; 150192bb6962SLisandro Dalcin row = ja + ia[i]; 150292bb6962SLisandro Dalcin for (j=0; j<len; j++) { 150392bb6962SLisandro Dalcin if (row[j] != i) { /* if not diagonal */ 150492bb6962SLisandro Dalcin jja[nnz+cnt++] = row[j]; 150592bb6962SLisandro Dalcin } 150692bb6962SLisandro Dalcin } 150792bb6962SLisandro Dalcin nnz += cnt; 150892bb6962SLisandro Dalcin iia[i+1] = nnz; 150992bb6962SLisandro Dalcin } 151092bb6962SLisandro Dalcin /* Partitioning of the adjacency matrix */ 15110298fd71SBarry Smith ierr = MatCreateMPIAdj(PETSC_COMM_SELF,na,na,iia,jja,NULL,&adj);CHKERRQ(ierr); 151292bb6962SLisandro Dalcin ierr = MatPartitioningSetAdjacency(mpart,adj);CHKERRQ(ierr); 151392bb6962SLisandro Dalcin ierr = MatPartitioningSetNParts(mpart,n);CHKERRQ(ierr); 151492bb6962SLisandro Dalcin ierr = MatPartitioningApply(mpart,&ispart);CHKERRQ(ierr); 151592bb6962SLisandro Dalcin ierr = ISPartitioningToNumbering(ispart,&isnumb);CHKERRQ(ierr); 1516fcfd50ebSBarry Smith ierr = MatDestroy(&adj);CHKERRQ(ierr); 151792bb6962SLisandro Dalcin foundpart = PETSC_TRUE; 151892bb6962SLisandro Dalcin } 151992bb6962SLisandro Dalcin ierr = MatRestoreRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr); 152092bb6962SLisandro Dalcin } 1521fcfd50ebSBarry Smith ierr = MatPartitioningDestroy(&mpart);CHKERRQ(ierr); 152292bb6962SLisandro Dalcin } 152392bb6962SLisandro Dalcin 1524785e854fSJed Brown ierr = PetscMalloc1(n,&is);CHKERRQ(ierr); 152592bb6962SLisandro Dalcin *outis = is; 152692bb6962SLisandro Dalcin 152792bb6962SLisandro Dalcin if (!foundpart) { 152892bb6962SLisandro Dalcin 152992bb6962SLisandro Dalcin /* Partitioning by contiguous chunks of rows */ 153092bb6962SLisandro Dalcin 153192bb6962SLisandro Dalcin PetscInt mbs = (rend-rstart)/bs; 153292bb6962SLisandro Dalcin PetscInt start = rstart; 153392bb6962SLisandro Dalcin for (i=0; i<n; i++) { 153492bb6962SLisandro Dalcin PetscInt count = (mbs/n + ((mbs % n) > i)) * bs; 153592bb6962SLisandro Dalcin ierr = ISCreateStride(PETSC_COMM_SELF,count,start,1,&is[i]);CHKERRQ(ierr); 153692bb6962SLisandro Dalcin start += count; 153792bb6962SLisandro Dalcin } 153892bb6962SLisandro Dalcin 153992bb6962SLisandro Dalcin } else { 154092bb6962SLisandro Dalcin 154192bb6962SLisandro Dalcin /* Partitioning by adjacency of diagonal block */ 154292bb6962SLisandro Dalcin 154392bb6962SLisandro Dalcin const PetscInt *numbering; 154492bb6962SLisandro Dalcin PetscInt *count,nidx,*indices,*newidx,start=0; 154592bb6962SLisandro Dalcin /* Get node count in each partition */ 1546785e854fSJed Brown ierr = PetscMalloc1(n,&count);CHKERRQ(ierr); 154792bb6962SLisandro Dalcin ierr = ISPartitioningCount(ispart,n,count);CHKERRQ(ierr); 154892bb6962SLisandro Dalcin if (isbaij && bs > 1) { /* adjust for the block-aij case */ 154992bb6962SLisandro Dalcin for (i=0; i<n; i++) count[i] *= bs; 155092bb6962SLisandro Dalcin } 155192bb6962SLisandro Dalcin /* Build indices from node numbering */ 155292bb6962SLisandro Dalcin ierr = ISGetLocalSize(isnumb,&nidx);CHKERRQ(ierr); 1553785e854fSJed Brown ierr = PetscMalloc1(nidx,&indices);CHKERRQ(ierr); 155492bb6962SLisandro Dalcin for (i=0; i<nidx; i++) indices[i] = i; /* needs to be initialized */ 155592bb6962SLisandro Dalcin ierr = ISGetIndices(isnumb,&numbering);CHKERRQ(ierr); 155692bb6962SLisandro Dalcin ierr = PetscSortIntWithPermutation(nidx,numbering,indices);CHKERRQ(ierr); 155792bb6962SLisandro Dalcin ierr = ISRestoreIndices(isnumb,&numbering);CHKERRQ(ierr); 155892bb6962SLisandro Dalcin if (isbaij && bs > 1) { /* adjust for the block-aij case */ 1559785e854fSJed Brown ierr = PetscMalloc1(nidx*bs,&newidx);CHKERRQ(ierr); 15602fa5cd67SKarl Rupp for (i=0; i<nidx; i++) { 15612fa5cd67SKarl Rupp for (j=0; j<bs; j++) newidx[i*bs+j] = indices[i]*bs + j; 15622fa5cd67SKarl Rupp } 156392bb6962SLisandro Dalcin ierr = PetscFree(indices);CHKERRQ(ierr); 156492bb6962SLisandro Dalcin nidx *= bs; 156592bb6962SLisandro Dalcin indices = newidx; 156692bb6962SLisandro Dalcin } 156792bb6962SLisandro Dalcin /* Shift to get global indices */ 156892bb6962SLisandro Dalcin for (i=0; i<nidx; i++) indices[i] += rstart; 156992bb6962SLisandro Dalcin 157092bb6962SLisandro Dalcin /* Build the index sets for each block */ 157192bb6962SLisandro Dalcin for (i=0; i<n; i++) { 157270b3c8c7SBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,count[i],&indices[start],PETSC_COPY_VALUES,&is[i]);CHKERRQ(ierr); 157392bb6962SLisandro Dalcin ierr = ISSort(is[i]);CHKERRQ(ierr); 157492bb6962SLisandro Dalcin start += count[i]; 157592bb6962SLisandro Dalcin } 157692bb6962SLisandro Dalcin 15773bf036e2SBarry Smith ierr = PetscFree(count);CHKERRQ(ierr); 15783bf036e2SBarry Smith ierr = PetscFree(indices);CHKERRQ(ierr); 1579fcfd50ebSBarry Smith ierr = ISDestroy(&isnumb);CHKERRQ(ierr); 1580fcfd50ebSBarry Smith ierr = ISDestroy(&ispart);CHKERRQ(ierr); 158192bb6962SLisandro Dalcin 158292bb6962SLisandro Dalcin } 158392bb6962SLisandro Dalcin PetscFunctionReturn(0); 158492bb6962SLisandro Dalcin } 158592bb6962SLisandro Dalcin 158692bb6962SLisandro Dalcin /*@C 158792bb6962SLisandro Dalcin PCASMDestroySubdomains - Destroys the index sets created with 158892bb6962SLisandro Dalcin PCASMCreateSubdomains(). Should be called after setting subdomains 158992bb6962SLisandro Dalcin with PCASMSetLocalSubdomains(). 159092bb6962SLisandro Dalcin 159192bb6962SLisandro Dalcin Collective 159292bb6962SLisandro Dalcin 159392bb6962SLisandro Dalcin Input Parameters: 159492bb6962SLisandro Dalcin + n - the number of index sets 15952b691e39SMatthew Knepley . is - the array of index sets 15960298fd71SBarry Smith - is_local - the array of local index sets, can be NULL 159792bb6962SLisandro Dalcin 159892bb6962SLisandro Dalcin Level: advanced 159992bb6962SLisandro Dalcin 160092bb6962SLisandro Dalcin .keywords: PC, ASM, additive Schwarz, create, subdomains, unstructured grid 160192bb6962SLisandro Dalcin 160292bb6962SLisandro Dalcin .seealso: PCASMCreateSubdomains(), PCASMSetLocalSubdomains() 160392bb6962SLisandro Dalcin @*/ 16047087cfbeSBarry Smith PetscErrorCode PCASMDestroySubdomains(PetscInt n, IS is[], IS is_local[]) 160592bb6962SLisandro Dalcin { 160692bb6962SLisandro Dalcin PetscInt i; 160792bb6962SLisandro Dalcin PetscErrorCode ierr; 16085fd66863SKarl Rupp 160992bb6962SLisandro Dalcin PetscFunctionBegin; 1610a35b7b57SDmitry Karpeev if (n <= 0) PetscFunctionReturn(0); 1611a35b7b57SDmitry Karpeev if (is) { 161292bb6962SLisandro Dalcin PetscValidPointer(is,2); 1613fcfd50ebSBarry Smith for (i=0; i<n; i++) { ierr = ISDestroy(&is[i]);CHKERRQ(ierr); } 161492bb6962SLisandro Dalcin ierr = PetscFree(is);CHKERRQ(ierr); 1615a35b7b57SDmitry Karpeev } 16162b691e39SMatthew Knepley if (is_local) { 16172b691e39SMatthew Knepley PetscValidPointer(is_local,3); 1618fcfd50ebSBarry Smith for (i=0; i<n; i++) { ierr = ISDestroy(&is_local[i]);CHKERRQ(ierr); } 16192b691e39SMatthew Knepley ierr = PetscFree(is_local);CHKERRQ(ierr); 16202b691e39SMatthew Knepley } 162192bb6962SLisandro Dalcin PetscFunctionReturn(0); 162292bb6962SLisandro Dalcin } 162392bb6962SLisandro Dalcin 16244b9ad928SBarry Smith /*@ 16254b9ad928SBarry Smith PCASMCreateSubdomains2D - Creates the index sets for the overlapping Schwarz 16264b9ad928SBarry Smith preconditioner for a two-dimensional problem on a regular grid. 16274b9ad928SBarry Smith 16284b9ad928SBarry Smith Not Collective 16294b9ad928SBarry Smith 16304b9ad928SBarry Smith Input Parameters: 16314b9ad928SBarry Smith + m, n - the number of mesh points in the x and y directions 16324b9ad928SBarry Smith . M, N - the number of subdomains in the x and y directions 16334b9ad928SBarry Smith . dof - degrees of freedom per node 16344b9ad928SBarry Smith - overlap - overlap in mesh lines 16354b9ad928SBarry Smith 16364b9ad928SBarry Smith Output Parameters: 16374b9ad928SBarry Smith + Nsub - the number of subdomains created 16383d03bb48SJed Brown . is - array of index sets defining overlapping (if overlap > 0) subdomains 16393d03bb48SJed Brown - is_local - array of index sets defining non-overlapping subdomains 16404b9ad928SBarry Smith 16414b9ad928SBarry Smith Note: 16424b9ad928SBarry Smith Presently PCAMSCreateSubdomains2d() is valid only for sequential 16434b9ad928SBarry Smith preconditioners. More general related routines are 16444b9ad928SBarry Smith PCASMSetTotalSubdomains() and PCASMSetLocalSubdomains(). 16454b9ad928SBarry Smith 16464b9ad928SBarry Smith Level: advanced 16474b9ad928SBarry Smith 16484b9ad928SBarry Smith .keywords: PC, ASM, additive Schwarz, create, subdomains, 2D, regular grid 16494b9ad928SBarry Smith 16504b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetLocalSubdomains(), PCASMGetSubKSP(), 16514b9ad928SBarry Smith PCASMSetOverlap() 16524b9ad928SBarry Smith @*/ 16537087cfbeSBarry Smith PetscErrorCode PCASMCreateSubdomains2D(PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt dof,PetscInt overlap,PetscInt *Nsub,IS **is,IS **is_local) 16544b9ad928SBarry Smith { 16553d03bb48SJed Brown PetscInt i,j,height,width,ystart,xstart,yleft,yright,xleft,xright,loc_outer; 16566849ba73SBarry Smith PetscErrorCode ierr; 165713f74950SBarry Smith PetscInt nidx,*idx,loc,ii,jj,count; 16584b9ad928SBarry Smith 16594b9ad928SBarry Smith PetscFunctionBegin; 1660e32f2f54SBarry Smith if (dof != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP," "); 16614b9ad928SBarry Smith 16624b9ad928SBarry Smith *Nsub = N*M; 1663854ce69bSBarry Smith ierr = PetscMalloc1(*Nsub,is);CHKERRQ(ierr); 1664854ce69bSBarry Smith ierr = PetscMalloc1(*Nsub,is_local);CHKERRQ(ierr); 16654b9ad928SBarry Smith ystart = 0; 16663d03bb48SJed Brown loc_outer = 0; 16674b9ad928SBarry Smith for (i=0; i<N; i++) { 16684b9ad928SBarry Smith height = n/N + ((n % N) > i); /* height of subdomain */ 1669e32f2f54SBarry Smith if (height < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many N subdomains for mesh dimension n"); 16704b9ad928SBarry Smith yleft = ystart - overlap; if (yleft < 0) yleft = 0; 16714b9ad928SBarry Smith yright = ystart + height + overlap; if (yright > n) yright = n; 16724b9ad928SBarry Smith xstart = 0; 16734b9ad928SBarry Smith for (j=0; j<M; j++) { 16744b9ad928SBarry Smith width = m/M + ((m % M) > j); /* width of subdomain */ 1675e32f2f54SBarry Smith if (width < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many M subdomains for mesh dimension m"); 16764b9ad928SBarry Smith xleft = xstart - overlap; if (xleft < 0) xleft = 0; 16774b9ad928SBarry Smith xright = xstart + width + overlap; if (xright > m) xright = m; 16784b9ad928SBarry Smith nidx = (xright - xleft)*(yright - yleft); 1679785e854fSJed Brown ierr = PetscMalloc1(nidx,&idx);CHKERRQ(ierr); 16804b9ad928SBarry Smith loc = 0; 16814b9ad928SBarry Smith for (ii=yleft; ii<yright; ii++) { 16824b9ad928SBarry Smith count = m*ii + xleft; 16832fa5cd67SKarl Rupp for (jj=xleft; jj<xright; jj++) idx[loc++] = count++; 16844b9ad928SBarry Smith } 168570b3c8c7SBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,nidx,idx,PETSC_COPY_VALUES,(*is)+loc_outer);CHKERRQ(ierr); 16863d03bb48SJed Brown if (overlap == 0) { 16873d03bb48SJed Brown ierr = PetscObjectReference((PetscObject)(*is)[loc_outer]);CHKERRQ(ierr); 16882fa5cd67SKarl Rupp 16893d03bb48SJed Brown (*is_local)[loc_outer] = (*is)[loc_outer]; 16903d03bb48SJed Brown } else { 16913d03bb48SJed Brown for (loc=0,ii=ystart; ii<ystart+height; ii++) { 16923d03bb48SJed Brown for (jj=xstart; jj<xstart+width; jj++) { 16933d03bb48SJed Brown idx[loc++] = m*ii + jj; 16943d03bb48SJed Brown } 16953d03bb48SJed Brown } 169670b3c8c7SBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,loc,idx,PETSC_COPY_VALUES,*is_local+loc_outer);CHKERRQ(ierr); 16973d03bb48SJed Brown } 16984b9ad928SBarry Smith ierr = PetscFree(idx);CHKERRQ(ierr); 16994b9ad928SBarry Smith xstart += width; 17003d03bb48SJed Brown loc_outer++; 17014b9ad928SBarry Smith } 17024b9ad928SBarry Smith ystart += height; 17034b9ad928SBarry Smith } 17044b9ad928SBarry Smith for (i=0; i<*Nsub; i++) { ierr = ISSort((*is)[i]);CHKERRQ(ierr); } 17054b9ad928SBarry Smith PetscFunctionReturn(0); 17064b9ad928SBarry Smith } 17074b9ad928SBarry Smith 17084b9ad928SBarry Smith /*@C 17094b9ad928SBarry Smith PCASMGetLocalSubdomains - Gets the local subdomains (for this processor 17104b9ad928SBarry Smith only) for the additive Schwarz preconditioner. 17114b9ad928SBarry Smith 1712ad4df100SBarry Smith Not Collective 17134b9ad928SBarry Smith 17144b9ad928SBarry Smith Input Parameter: 17154b9ad928SBarry Smith . pc - the preconditioner context 17164b9ad928SBarry Smith 17174b9ad928SBarry Smith Output Parameters: 17184b9ad928SBarry Smith + n - the number of subdomains for this processor (default value = 1) 17192b691e39SMatthew Knepley . is - the index sets that define the subdomains for this processor 17200298fd71SBarry Smith - is_local - the index sets that define the local part of the subdomains for this processor (can be NULL) 17214b9ad928SBarry Smith 17224b9ad928SBarry Smith 17234b9ad928SBarry Smith Notes: 17244b9ad928SBarry Smith The IS numbering is in the parallel, global numbering of the vector. 17254b9ad928SBarry Smith 17264b9ad928SBarry Smith Level: advanced 17274b9ad928SBarry Smith 17284b9ad928SBarry Smith .keywords: PC, ASM, set, local, subdomains, additive Schwarz 17294b9ad928SBarry Smith 17304b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(), 17314b9ad928SBarry Smith PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubmatrices() 17324b9ad928SBarry Smith @*/ 17337087cfbeSBarry Smith PetscErrorCode PCASMGetLocalSubdomains(PC pc,PetscInt *n,IS *is[],IS *is_local[]) 17344b9ad928SBarry Smith { 17352a808120SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 173692bb6962SLisandro Dalcin PetscErrorCode ierr; 1737ace3abfcSBarry Smith PetscBool match; 17384b9ad928SBarry Smith 17394b9ad928SBarry Smith PetscFunctionBegin; 17400700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 17414482741eSBarry Smith PetscValidIntPointer(n,2); 174292bb6962SLisandro Dalcin if (is) PetscValidPointer(is,3); 1743251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)pc,PCASM,&match);CHKERRQ(ierr); 17442a808120SBarry Smith if (!match) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"PC is not a PCASM"); 17454b9ad928SBarry Smith if (n) *n = osm->n_local_true; 17464b9ad928SBarry Smith if (is) *is = osm->is; 17472b691e39SMatthew Knepley if (is_local) *is_local = osm->is_local; 17484b9ad928SBarry Smith PetscFunctionReturn(0); 17494b9ad928SBarry Smith } 17504b9ad928SBarry Smith 17514b9ad928SBarry Smith /*@C 17524b9ad928SBarry Smith PCASMGetLocalSubmatrices - Gets the local submatrices (for this processor 17534b9ad928SBarry Smith only) for the additive Schwarz preconditioner. 17544b9ad928SBarry Smith 1755ad4df100SBarry Smith Not Collective 17564b9ad928SBarry Smith 17574b9ad928SBarry Smith Input Parameter: 17584b9ad928SBarry Smith . pc - the preconditioner context 17594b9ad928SBarry Smith 17604b9ad928SBarry Smith Output Parameters: 17614b9ad928SBarry Smith + n - the number of matrices for this processor (default value = 1) 17624b9ad928SBarry Smith - mat - the matrices 17634b9ad928SBarry Smith 17644b9ad928SBarry Smith Level: advanced 17654b9ad928SBarry Smith 1766cf739d55SBarry Smith Notes: Call after PCSetUp() (or KSPSetUp()) but before PCApply() (or KSPApply()) and before PCSetUpOnBlocks()) 1767cf739d55SBarry Smith 1768cf739d55SBarry Smith Usually one would use PCSetModifySubmatrices() to change the submatrices in building the preconditioner. 1769cf739d55SBarry Smith 17704b9ad928SBarry Smith .keywords: PC, ASM, set, local, subdomains, additive Schwarz, block Jacobi 17714b9ad928SBarry Smith 17724b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(), 1773cf739d55SBarry Smith PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubdomains(), PCSetModifySubmatrices() 17744b9ad928SBarry Smith @*/ 17757087cfbeSBarry Smith PetscErrorCode PCASMGetLocalSubmatrices(PC pc,PetscInt *n,Mat *mat[]) 17764b9ad928SBarry Smith { 17774b9ad928SBarry Smith PC_ASM *osm; 177892bb6962SLisandro Dalcin PetscErrorCode ierr; 1779ace3abfcSBarry Smith PetscBool match; 17804b9ad928SBarry Smith 17814b9ad928SBarry Smith PetscFunctionBegin; 17820700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 178392bb6962SLisandro Dalcin PetscValidIntPointer(n,2); 178492bb6962SLisandro Dalcin if (mat) PetscValidPointer(mat,3); 1785ce94432eSBarry Smith if (!pc->setupcalled) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call after KSPSetUP() or PCSetUp()."); 1786251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)pc,PCASM,&match);CHKERRQ(ierr); 178792bb6962SLisandro Dalcin if (!match) { 178892bb6962SLisandro Dalcin if (n) *n = 0; 17890298fd71SBarry Smith if (mat) *mat = NULL; 179092bb6962SLisandro Dalcin } else { 17914b9ad928SBarry Smith osm = (PC_ASM*)pc->data; 17924b9ad928SBarry Smith if (n) *n = osm->n_local_true; 17934b9ad928SBarry Smith if (mat) *mat = osm->pmat; 179492bb6962SLisandro Dalcin } 17954b9ad928SBarry Smith PetscFunctionReturn(0); 17964b9ad928SBarry Smith } 1797d709ea83SDmitry Karpeev 1798d709ea83SDmitry Karpeev /*@ 1799d709ea83SDmitry Karpeev PCASMSetDMSubdomains - Indicates whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible. 1800f1ee410cSBarry Smith 1801d709ea83SDmitry Karpeev Logically Collective 1802d709ea83SDmitry Karpeev 1803d709ea83SDmitry Karpeev Input Parameter: 1804d709ea83SDmitry Karpeev + pc - the preconditioner 1805d709ea83SDmitry Karpeev - flg - boolean indicating whether to use subdomains defined by the DM 1806d709ea83SDmitry Karpeev 1807d709ea83SDmitry Karpeev Options Database Key: 1808d709ea83SDmitry Karpeev . -pc_asm_dm_subdomains 1809d709ea83SDmitry Karpeev 1810d709ea83SDmitry Karpeev Level: intermediate 1811d709ea83SDmitry Karpeev 1812d709ea83SDmitry Karpeev Notes: 1813d709ea83SDmitry Karpeev PCASMSetTotalSubdomains() and PCASMSetOverlap() take precedence over PCASMSetDMSubdomains(), 1814d709ea83SDmitry Karpeev so setting either of the first two effectively turns the latter off. 1815d709ea83SDmitry Karpeev 1816d709ea83SDmitry Karpeev .keywords: PC, ASM, DM, set, subdomains, additive Schwarz 1817d709ea83SDmitry Karpeev 1818d709ea83SDmitry Karpeev .seealso: PCASMGetDMSubdomains(), PCASMSetTotalSubdomains(), PCASMSetOverlap() 1819d709ea83SDmitry Karpeev PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubdomains() 1820d709ea83SDmitry Karpeev @*/ 1821d709ea83SDmitry Karpeev PetscErrorCode PCASMSetDMSubdomains(PC pc,PetscBool flg) 1822d709ea83SDmitry Karpeev { 1823d709ea83SDmitry Karpeev PC_ASM *osm = (PC_ASM*)pc->data; 1824d709ea83SDmitry Karpeev PetscErrorCode ierr; 1825d709ea83SDmitry Karpeev PetscBool match; 1826d709ea83SDmitry Karpeev 1827d709ea83SDmitry Karpeev PetscFunctionBegin; 1828d709ea83SDmitry Karpeev PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1829d709ea83SDmitry Karpeev PetscValidLogicalCollectiveBool(pc,flg,2); 1830d709ea83SDmitry Karpeev if (pc->setupcalled) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Not for a setup PC."); 1831d709ea83SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)pc,PCASM,&match);CHKERRQ(ierr); 1832d709ea83SDmitry Karpeev if (match) { 1833d709ea83SDmitry Karpeev osm->dm_subdomains = flg; 1834d709ea83SDmitry Karpeev } 1835d709ea83SDmitry Karpeev PetscFunctionReturn(0); 1836d709ea83SDmitry Karpeev } 1837d709ea83SDmitry Karpeev 1838d709ea83SDmitry Karpeev /*@ 1839d709ea83SDmitry Karpeev PCASMGetDMSubdomains - Returns flag indicating whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible. 1840d709ea83SDmitry Karpeev Not Collective 1841d709ea83SDmitry Karpeev 1842d709ea83SDmitry Karpeev Input Parameter: 1843d709ea83SDmitry Karpeev . pc - the preconditioner 1844d709ea83SDmitry Karpeev 1845d709ea83SDmitry Karpeev Output Parameter: 1846d709ea83SDmitry Karpeev . flg - boolean indicating whether to use subdomains defined by the DM 1847d709ea83SDmitry Karpeev 1848d709ea83SDmitry Karpeev Level: intermediate 1849d709ea83SDmitry Karpeev 1850d709ea83SDmitry Karpeev .keywords: PC, ASM, DM, set, subdomains, additive Schwarz 1851d709ea83SDmitry Karpeev 1852d709ea83SDmitry Karpeev .seealso: PCASMSetDMSubdomains(), PCASMSetTotalSubdomains(), PCASMSetOverlap() 1853d709ea83SDmitry Karpeev PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubdomains() 1854d709ea83SDmitry Karpeev @*/ 1855d709ea83SDmitry Karpeev PetscErrorCode PCASMGetDMSubdomains(PC pc,PetscBool* flg) 1856d709ea83SDmitry Karpeev { 1857d709ea83SDmitry Karpeev PC_ASM *osm = (PC_ASM*)pc->data; 1858d709ea83SDmitry Karpeev PetscErrorCode ierr; 1859d709ea83SDmitry Karpeev PetscBool match; 1860d709ea83SDmitry Karpeev 1861d709ea83SDmitry Karpeev PetscFunctionBegin; 1862d709ea83SDmitry Karpeev PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1863d709ea83SDmitry Karpeev PetscValidPointer(flg,2); 1864d709ea83SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)pc,PCASM,&match);CHKERRQ(ierr); 1865d709ea83SDmitry Karpeev if (match) { 1866d709ea83SDmitry Karpeev if (flg) *flg = osm->dm_subdomains; 1867d709ea83SDmitry Karpeev } 1868d709ea83SDmitry Karpeev PetscFunctionReturn(0); 1869d709ea83SDmitry Karpeev } 187080ec0b7dSPatrick Sanan 187180ec0b7dSPatrick Sanan /*@ 187280ec0b7dSPatrick Sanan PCASMGetSubMatType - Gets the matrix type used for ASM subsolves, as a string. 187380ec0b7dSPatrick Sanan 187480ec0b7dSPatrick Sanan Not Collective 187580ec0b7dSPatrick Sanan 187680ec0b7dSPatrick Sanan Input Parameter: 187780ec0b7dSPatrick Sanan . pc - the PC 187880ec0b7dSPatrick Sanan 187980ec0b7dSPatrick Sanan Output Parameter: 1880f1ee410cSBarry Smith . -pc_asm_sub_mat_type - name of matrix type 188180ec0b7dSPatrick Sanan 188280ec0b7dSPatrick Sanan Level: advanced 188380ec0b7dSPatrick Sanan 188480ec0b7dSPatrick Sanan .keywords: PC, PCASM, MatType, set 188580ec0b7dSPatrick Sanan 188680ec0b7dSPatrick Sanan .seealso: PCASMSetSubMatType(), PCASM, PCSetType(), VecSetType(), MatType, Mat 188780ec0b7dSPatrick Sanan @*/ 188880ec0b7dSPatrick Sanan PetscErrorCode PCASMGetSubMatType(PC pc,MatType *sub_mat_type){ 188980ec0b7dSPatrick Sanan PetscErrorCode ierr; 189080ec0b7dSPatrick Sanan 189180ec0b7dSPatrick Sanan ierr = PetscTryMethod(pc,"PCASMGetSubMatType_C",(PC,MatType*),(pc,sub_mat_type));CHKERRQ(ierr); 189280ec0b7dSPatrick Sanan PetscFunctionReturn(0); 189380ec0b7dSPatrick Sanan } 189480ec0b7dSPatrick Sanan 189580ec0b7dSPatrick Sanan /*@ 189680ec0b7dSPatrick Sanan PCASMSetSubMatType - Set the type of matrix used for ASM subsolves 189780ec0b7dSPatrick Sanan 189880ec0b7dSPatrick Sanan Collective on Mat 189980ec0b7dSPatrick Sanan 190080ec0b7dSPatrick Sanan Input Parameters: 190180ec0b7dSPatrick Sanan + pc - the PC object 190280ec0b7dSPatrick Sanan - sub_mat_type - matrix type 190380ec0b7dSPatrick Sanan 190480ec0b7dSPatrick Sanan Options Database Key: 190580ec0b7dSPatrick 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. 190680ec0b7dSPatrick Sanan 190780ec0b7dSPatrick Sanan Notes: 190880ec0b7dSPatrick Sanan See "${PETSC_DIR}/include/petscmat.h" for available types 190980ec0b7dSPatrick Sanan 191080ec0b7dSPatrick Sanan Level: advanced 191180ec0b7dSPatrick Sanan 191280ec0b7dSPatrick Sanan .keywords: PC, PCASM, MatType, set 191380ec0b7dSPatrick Sanan 191480ec0b7dSPatrick Sanan .seealso: PCASMGetSubMatType(), PCASM, PCSetType(), VecSetType(), MatType, Mat 191580ec0b7dSPatrick Sanan @*/ 191680ec0b7dSPatrick Sanan PetscErrorCode PCASMSetSubMatType(PC pc,MatType sub_mat_type) 191780ec0b7dSPatrick Sanan { 191880ec0b7dSPatrick Sanan PetscErrorCode ierr; 191980ec0b7dSPatrick Sanan 192080ec0b7dSPatrick Sanan ierr = PetscTryMethod(pc,"PCASMSetSubMatType_C",(PC,MatType),(pc,sub_mat_type));CHKERRQ(ierr); 192180ec0b7dSPatrick Sanan PetscFunctionReturn(0); 192280ec0b7dSPatrick Sanan } 1923