14b9ad928SBarry Smith /* 24b9ad928SBarry Smith This file defines an additive Schwarz preconditioner for any Mat implementation. 34b9ad928SBarry Smith 44b9ad928SBarry Smith Note that each processor may have any number of subdomains. But in order to 54b9ad928SBarry Smith deal easily with the VecScatter(), we treat each processor as if it has the 64b9ad928SBarry Smith same number of subdomains. 74b9ad928SBarry Smith 84b9ad928SBarry Smith n - total number of true subdomains on all processors 94b9ad928SBarry Smith n_local_true - actual number of subdomains on this processor 104b9ad928SBarry Smith n_local = maximum over all processors of n_local_true 114b9ad928SBarry Smith */ 12af0996ceSBarry Smith #include <petsc/private/pcimpl.h> /*I "petscpc.h" I*/ 131e25c274SJed Brown #include <petscdm.h> 144b9ad928SBarry Smith 154b9ad928SBarry Smith typedef struct { 1613f74950SBarry Smith PetscInt n, n_local, n_local_true; 1713f74950SBarry Smith PetscInt overlap; /* overlap requested by user */ 184b9ad928SBarry Smith KSP *ksp; /* linear solvers for each block */ 19b0de9f23SBarry Smith VecScatter restriction; /* mapping from global to overlapping (process) subdomain*/ 20b0de9f23SBarry Smith VecScatter *lrestriction; /* mapping from subregion to overlapping (process) subdomain */ 21b0de9f23SBarry Smith VecScatter *lprolongation; /* mapping from non-overlapping subregion to overlapping (process) subdomain; used for restrict additive version of algorithms */ 221dd8081eSeaulisa Vec lx, ly; /* work vectors */ 231dd8081eSeaulisa Vec *x,*y; /* work vectors */ 241dd8081eSeaulisa IS lis; /* index set that defines each overlapping multiplicative (process) subdomain */ 252b691e39SMatthew Knepley IS *is; /* index set that defines each overlapping subdomain */ 262b691e39SMatthew Knepley IS *is_local; /* index set that defines each non-overlapping subdomain, may be NULL */ 274b9ad928SBarry Smith Mat *mat,*pmat; /* mat is not currently used */ 284b9ad928SBarry Smith PCASMType type; /* use reduced interpolation, restriction or both */ 29ace3abfcSBarry Smith PetscBool type_set; /* if user set this value (so won't change it for symmetric problems) */ 30ace3abfcSBarry Smith PetscBool same_local_solves; /* flag indicating whether all local solvers are same */ 31ace3abfcSBarry Smith PetscBool sort_indices; /* flag to sort subdomain indices */ 32d709ea83SDmitry Karpeev PetscBool dm_subdomains; /* whether DM is allowed to define subdomains */ 3312cd4985SMatthew G. Knepley PCCompositeType loctype; /* the type of composition for local solves */ 3480ec0b7dSPatrick Sanan MatType sub_mat_type; /* the type of Mat used for subdomain solves (can be MATSAME or NULL) */ 35fb745f2cSMatthew G. Knepley /* For multiplicative solve */ 36235cc4ceSMatthew G. Knepley Mat *lmats; /* submatrices for overlapping multiplicative (process) subdomain */ 374b9ad928SBarry Smith } PC_ASM; 384b9ad928SBarry Smith 396849ba73SBarry Smith static PetscErrorCode PCView_ASM(PC pc,PetscViewer viewer) 404b9ad928SBarry Smith { 4192bb6962SLisandro Dalcin PC_ASM *osm = (PC_ASM*)pc->data; 426849ba73SBarry Smith PetscErrorCode ierr; 4313f74950SBarry Smith PetscMPIInt rank; 444d219a6aSLisandro Dalcin PetscInt i,bsz; 45ace3abfcSBarry Smith PetscBool iascii,isstring; 464b9ad928SBarry Smith PetscViewer sviewer; 474b9ad928SBarry Smith 484b9ad928SBarry Smith PetscFunctionBegin; 49251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 50251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr); 5132077d6dSBarry Smith if (iascii) { 523d03bb48SJed Brown char overlaps[256] = "user-defined overlap",blocks[256] = "total subdomain blocks not yet set"; 538caf3d72SBarry Smith if (osm->overlap >= 0) {ierr = PetscSNPrintf(overlaps,sizeof(overlaps),"amount of overlap = %D",osm->overlap);CHKERRQ(ierr);} 548caf3d72SBarry Smith if (osm->n > 0) {ierr = PetscSNPrintf(blocks,sizeof(blocks),"total subdomain blocks = %D",osm->n);CHKERRQ(ierr);} 55efd4aadfSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %s, %s\n",blocks,overlaps);CHKERRQ(ierr); 56efd4aadfSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," restriction/interpolation type - %s\n",PCASMTypes[osm->type]);CHKERRQ(ierr); 57e64f0791SPatrick Sanan if (osm->dm_subdomains) {ierr = PetscViewerASCIIPrintf(viewer," Additive Schwarz: using DM to define subdomains\n");CHKERRQ(ierr);} 5812cd4985SMatthew G. Knepley if (osm->loctype != PC_COMPOSITE_ADDITIVE) {ierr = PetscViewerASCIIPrintf(viewer," Additive Schwarz: local solve composition type - %s\n",PCCompositeTypes[osm->loctype]);CHKERRQ(ierr);} 59ce94432eSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);CHKERRQ(ierr); 6092bb6962SLisandro Dalcin if (osm->same_local_solves) { 614d219a6aSLisandro Dalcin if (osm->ksp) { 624b9ad928SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Local solve is same for all blocks, in the following KSP and PC objects:\n");CHKERRQ(ierr); 63c5e4d11fSDmitry Karpeev ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 644d219a6aSLisandro Dalcin if (!rank) { 654b9ad928SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 6692bb6962SLisandro Dalcin ierr = KSPView(osm->ksp[0],sviewer);CHKERRQ(ierr); 674b9ad928SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 684b9ad928SBarry Smith } 69c5e4d11fSDmitry Karpeev ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 704d219a6aSLisandro Dalcin } 714b9ad928SBarry Smith } else { 72c5e4d11fSDmitry Karpeev ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 734d219a6aSLisandro Dalcin ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] number of local blocks = %D\n",(int)rank,osm->n_local_true);CHKERRQ(ierr); 744d219a6aSLisandro Dalcin ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 754b9ad928SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Local solve info for each block is in the following KSP and PC objects:\n");CHKERRQ(ierr); 764b9ad928SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 774d219a6aSLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr); 78c5e4d11fSDmitry Karpeev ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 79dd2fa690SBarry Smith for (i=0; i<osm->n_local_true; i++) { 804d219a6aSLisandro Dalcin ierr = ISGetLocalSize(osm->is[i],&bsz);CHKERRQ(ierr); 814d219a6aSLisandro Dalcin ierr = PetscViewerASCIISynchronizedPrintf(sviewer,"[%d] local block number %D, size = %D\n",(int)rank,i,bsz);CHKERRQ(ierr); 8292bb6962SLisandro Dalcin ierr = KSPView(osm->ksp[i],sviewer);CHKERRQ(ierr); 834d219a6aSLisandro Dalcin ierr = PetscViewerASCIISynchronizedPrintf(sviewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr); 844b9ad928SBarry Smith } 85c5e4d11fSDmitry Karpeev ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 864b9ad928SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 874b9ad928SBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 88c5e4d11fSDmitry Karpeev ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 894b9ad928SBarry Smith } 904b9ad928SBarry Smith } else if (isstring) { 914d219a6aSLisandro Dalcin ierr = PetscViewerStringSPrintf(viewer," blocks=%D, overlap=%D, type=%s",osm->n,osm->overlap,PCASMTypes[osm->type]);CHKERRQ(ierr); 92c5e4d11fSDmitry Karpeev ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 9392bb6962SLisandro Dalcin if (osm->ksp) {ierr = KSPView(osm->ksp[0],sviewer);CHKERRQ(ierr);} 94c5e4d11fSDmitry Karpeev ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 954b9ad928SBarry Smith } 964b9ad928SBarry Smith PetscFunctionReturn(0); 974b9ad928SBarry Smith } 984b9ad928SBarry Smith 9992bb6962SLisandro Dalcin static PetscErrorCode PCASMPrintSubdomains(PC pc) 10092bb6962SLisandro Dalcin { 10192bb6962SLisandro Dalcin PC_ASM *osm = (PC_ASM*)pc->data; 10292bb6962SLisandro Dalcin const char *prefix; 10392bb6962SLisandro Dalcin char fname[PETSC_MAX_PATH_LEN+1]; 104643535aeSDmitry Karpeev PetscViewer viewer, sviewer; 105643535aeSDmitry Karpeev char *s; 10692bb6962SLisandro Dalcin PetscInt i,j,nidx; 10792bb6962SLisandro Dalcin const PetscInt *idx; 108643535aeSDmitry Karpeev PetscMPIInt rank, size; 10992bb6962SLisandro Dalcin PetscErrorCode ierr; 110846783a0SBarry Smith 11192bb6962SLisandro Dalcin PetscFunctionBegin; 112ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size);CHKERRQ(ierr); 113ce94432eSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank);CHKERRQ(ierr); 11492bb6962SLisandro Dalcin ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 115c5929fdfSBarry Smith ierr = PetscOptionsGetString(NULL,prefix,"-pc_asm_print_subdomains",fname,PETSC_MAX_PATH_LEN,NULL);CHKERRQ(ierr); 11692bb6962SLisandro Dalcin if (fname[0] == 0) { ierr = PetscStrcpy(fname,"stdout");CHKERRQ(ierr); }; 117ce94432eSBarry Smith ierr = PetscViewerASCIIOpen(PetscObjectComm((PetscObject)pc),fname,&viewer);CHKERRQ(ierr); 118643535aeSDmitry Karpeev for (i=0; i<osm->n_local; i++) { 119643535aeSDmitry Karpeev if (i < osm->n_local_true) { 12092bb6962SLisandro Dalcin ierr = ISGetLocalSize(osm->is[i],&nidx);CHKERRQ(ierr); 12192bb6962SLisandro Dalcin ierr = ISGetIndices(osm->is[i],&idx);CHKERRQ(ierr); 122643535aeSDmitry Karpeev /* Print to a string viewer; no more than 15 characters per index plus 512 char for the header.*/ 123854ce69bSBarry Smith ierr = PetscMalloc1(16*(nidx+1)+512, &s);CHKERRQ(ierr); 124643535aeSDmitry Karpeev ierr = PetscViewerStringOpen(PETSC_COMM_SELF, s, 16*(nidx+1)+512, &sviewer);CHKERRQ(ierr); 125643535aeSDmitry Karpeev ierr = PetscViewerStringSPrintf(sviewer, "[%D:%D] Subdomain %D with overlap:\n", rank, size, i);CHKERRQ(ierr); 12692bb6962SLisandro Dalcin for (j=0; j<nidx; j++) { 127643535aeSDmitry Karpeev ierr = PetscViewerStringSPrintf(sviewer,"%D ",idx[j]);CHKERRQ(ierr); 12892bb6962SLisandro Dalcin } 12992bb6962SLisandro Dalcin ierr = ISRestoreIndices(osm->is[i],&idx);CHKERRQ(ierr); 130643535aeSDmitry Karpeev ierr = PetscViewerStringSPrintf(sviewer,"\n");CHKERRQ(ierr); 131643535aeSDmitry Karpeev ierr = PetscViewerDestroy(&sviewer);CHKERRQ(ierr); 132c5e4d11fSDmitry Karpeev ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 133643535aeSDmitry Karpeev ierr = PetscViewerASCIISynchronizedPrintf(viewer, s);CHKERRQ(ierr); 134643535aeSDmitry Karpeev ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 135c5e4d11fSDmitry Karpeev ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 136643535aeSDmitry Karpeev ierr = PetscFree(s);CHKERRQ(ierr); 1372b691e39SMatthew Knepley if (osm->is_local) { 138643535aeSDmitry Karpeev /* Print to a string viewer; no more than 15 characters per index plus 512 char for the header.*/ 139854ce69bSBarry Smith ierr = PetscMalloc1(16*(nidx+1)+512, &s);CHKERRQ(ierr); 140643535aeSDmitry Karpeev ierr = PetscViewerStringOpen(PETSC_COMM_SELF, s, 16*(nidx+1)+512, &sviewer);CHKERRQ(ierr); 14109d011a0SDmitry Karpeev ierr = PetscViewerStringSPrintf(sviewer, "[%D:%D] Subdomain %D without overlap:\n", rank, size, i);CHKERRQ(ierr); 1422b691e39SMatthew Knepley ierr = ISGetLocalSize(osm->is_local[i],&nidx);CHKERRQ(ierr); 1432b691e39SMatthew Knepley ierr = ISGetIndices(osm->is_local[i],&idx);CHKERRQ(ierr); 1442b691e39SMatthew Knepley for (j=0; j<nidx; j++) { 14509d011a0SDmitry Karpeev ierr = PetscViewerStringSPrintf(sviewer,"%D ",idx[j]);CHKERRQ(ierr); 1462b691e39SMatthew Knepley } 1472b691e39SMatthew Knepley ierr = ISRestoreIndices(osm->is_local[i],&idx);CHKERRQ(ierr); 14809d011a0SDmitry Karpeev ierr = PetscViewerStringSPrintf(sviewer,"\n");CHKERRQ(ierr); 149643535aeSDmitry Karpeev ierr = PetscViewerDestroy(&sviewer);CHKERRQ(ierr); 150c5e4d11fSDmitry Karpeev ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 151643535aeSDmitry Karpeev ierr = PetscViewerASCIISynchronizedPrintf(viewer, s);CHKERRQ(ierr); 152643535aeSDmitry Karpeev ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 153c5e4d11fSDmitry Karpeev ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 154643535aeSDmitry Karpeev ierr = PetscFree(s);CHKERRQ(ierr); 155643535aeSDmitry Karpeev } 1562fa5cd67SKarl Rupp } else { 157643535aeSDmitry Karpeev /* Participate in collective viewer calls. */ 158c5e4d11fSDmitry Karpeev ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 159643535aeSDmitry Karpeev ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 160c5e4d11fSDmitry Karpeev ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 161643535aeSDmitry Karpeev /* Assume either all ranks have is_local or none do. */ 162643535aeSDmitry Karpeev if (osm->is_local) { 163c5e4d11fSDmitry Karpeev ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 164643535aeSDmitry Karpeev ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 165c5e4d11fSDmitry Karpeev ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 166643535aeSDmitry Karpeev } 167fdc48646SDmitry Karpeev } 16892bb6962SLisandro Dalcin } 16992bb6962SLisandro Dalcin ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 170fcfd50ebSBarry Smith ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 17192bb6962SLisandro Dalcin PetscFunctionReturn(0); 17292bb6962SLisandro Dalcin } 17392bb6962SLisandro Dalcin 1746849ba73SBarry Smith static PetscErrorCode PCSetUp_ASM(PC pc) 1754b9ad928SBarry Smith { 1764b9ad928SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 1776849ba73SBarry Smith PetscErrorCode ierr; 178ace3abfcSBarry Smith PetscBool symset,flg; 17987e86712SBarry Smith PetscInt i,m,m_local; 1804b9ad928SBarry Smith MatReuse scall = MAT_REUSE_MATRIX; 1814b9ad928SBarry Smith IS isl; 1824b9ad928SBarry Smith KSP ksp; 1834b9ad928SBarry Smith PC subpc; 1842dcb1b2aSMatthew Knepley const char *prefix,*pprefix; 18523ce1328SBarry Smith Vec vec; 1860298fd71SBarry Smith DM *domain_dm = NULL; 1874b9ad928SBarry Smith 1884b9ad928SBarry Smith PetscFunctionBegin; 1894b9ad928SBarry Smith if (!pc->setupcalled) { 190265a2120SBarry Smith PetscInt m; 19192bb6962SLisandro Dalcin 19292bb6962SLisandro Dalcin if (!osm->type_set) { 19392bb6962SLisandro Dalcin ierr = MatIsSymmetricKnown(pc->pmat,&symset,&flg);CHKERRQ(ierr); 1942fa5cd67SKarl Rupp if (symset && flg) osm->type = PC_ASM_BASIC; 19592bb6962SLisandro Dalcin } 19692bb6962SLisandro Dalcin 1972b837212SDmitry Karpeev /* Note: if subdomains have been set either via PCASMSetTotalSubdomains() or via PCASMSetLocalSubdomains(), osm->n_local_true will not be PETSC_DECIDE */ 1982b837212SDmitry Karpeev if (osm->n_local_true == PETSC_DECIDE) { 19969ca1f37SDmitry Karpeev /* no subdomains given */ 20065db9045SDmitry Karpeev /* try pc->dm first, if allowed */ 201d709ea83SDmitry Karpeev if (osm->dm_subdomains && pc->dm) { 202feb221f8SDmitry Karpeev PetscInt num_domains, d; 203feb221f8SDmitry Karpeev char **domain_names; 2048d4ac253SDmitry Karpeev IS *inner_domain_is, *outer_domain_is; 2058d4ac253SDmitry Karpeev ierr = DMCreateDomainDecomposition(pc->dm, &num_domains, &domain_names, &inner_domain_is, &outer_domain_is, &domain_dm);CHKERRQ(ierr); 206704f0395SPatrick Sanan osm->overlap = -1; /* We do not want to increase the overlap of the IS. 207704f0395SPatrick Sanan A future improvement of this code might allow one to use 208704f0395SPatrick Sanan DM-defined subdomains and also increase the overlap, 209704f0395SPatrick Sanan but that is not currently supported */ 21069ca1f37SDmitry Karpeev if (num_domains) { 2118d4ac253SDmitry Karpeev ierr = PCASMSetLocalSubdomains(pc, num_domains, outer_domain_is, inner_domain_is);CHKERRQ(ierr); 21269ca1f37SDmitry Karpeev } 213feb221f8SDmitry Karpeev for (d = 0; d < num_domains; ++d) { 214a35b7b57SDmitry Karpeev if (domain_names) {ierr = PetscFree(domain_names[d]);CHKERRQ(ierr);} 215a35b7b57SDmitry Karpeev if (inner_domain_is) {ierr = ISDestroy(&inner_domain_is[d]);CHKERRQ(ierr);} 216a35b7b57SDmitry Karpeev if (outer_domain_is) {ierr = ISDestroy(&outer_domain_is[d]);CHKERRQ(ierr);} 217feb221f8SDmitry Karpeev } 218feb221f8SDmitry Karpeev ierr = PetscFree(domain_names);CHKERRQ(ierr); 2198d4ac253SDmitry Karpeev ierr = PetscFree(inner_domain_is);CHKERRQ(ierr); 2208d4ac253SDmitry Karpeev ierr = PetscFree(outer_domain_is);CHKERRQ(ierr); 221feb221f8SDmitry Karpeev } 2222b837212SDmitry Karpeev if (osm->n_local_true == PETSC_DECIDE) { 22369ca1f37SDmitry Karpeev /* still no subdomains; use one subdomain per processor */ 2242b837212SDmitry Karpeev osm->n_local_true = 1; 225feb221f8SDmitry Karpeev } 2262b837212SDmitry Karpeev } 2272b837212SDmitry Karpeev { /* determine the global and max number of subdomains */ 2286ac3741eSJed Brown struct {PetscInt max,sum;} inwork,outwork; 229c10200c1SHong Zhang PetscMPIInt size; 230c10200c1SHong Zhang 2316ac3741eSJed Brown inwork.max = osm->n_local_true; 2326ac3741eSJed Brown inwork.sum = osm->n_local_true; 233367daffbSBarry Smith ierr = MPIU_Allreduce(&inwork,&outwork,1,MPIU_2INT,MPIU_MAXSUM_OP,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr); 2346ac3741eSJed Brown osm->n_local = outwork.max; 2356ac3741eSJed Brown osm->n = outwork.sum; 236c10200c1SHong Zhang 237c10200c1SHong Zhang ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr); 238c10200c1SHong Zhang if (outwork.max == 1 && outwork.sum == size) { 2397dae84e0SHong Zhang /* osm->n_local_true = 1 on all processes, set this option may enable use of optimized MatCreateSubMatrices() implementation */ 240c10200c1SHong Zhang ierr = MatSetOption(pc->pmat,MAT_SUBMAT_SINGLEIS,PETSC_TRUE);CHKERRQ(ierr); 241c10200c1SHong Zhang } 2424b9ad928SBarry Smith } 24378904715SLisandro Dalcin if (!osm->is) { /* create the index sets */ 24478904715SLisandro Dalcin ierr = PCASMCreateSubdomains(pc->pmat,osm->n_local_true,&osm->is);CHKERRQ(ierr); 2454b9ad928SBarry Smith } 246f5234e35SJed Brown if (osm->n_local_true > 1 && !osm->is_local) { 247785e854fSJed Brown ierr = PetscMalloc1(osm->n_local_true,&osm->is_local);CHKERRQ(ierr); 248f5234e35SJed Brown for (i=0; i<osm->n_local_true; i++) { 249f5234e35SJed Brown if (osm->overlap > 0) { /* With positive overlap, osm->is[i] will be modified */ 250f5234e35SJed Brown ierr = ISDuplicate(osm->is[i],&osm->is_local[i]);CHKERRQ(ierr); 251f5234e35SJed Brown ierr = ISCopy(osm->is[i],osm->is_local[i]);CHKERRQ(ierr); 252f5234e35SJed Brown } else { 253f5234e35SJed Brown ierr = PetscObjectReference((PetscObject)osm->is[i]);CHKERRQ(ierr); 254f5234e35SJed Brown osm->is_local[i] = osm->is[i]; 255f5234e35SJed Brown } 256f5234e35SJed Brown } 257f5234e35SJed Brown } 25892bb6962SLisandro Dalcin ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 25990d69ab7SBarry Smith flg = PETSC_FALSE; 260c5929fdfSBarry Smith ierr = PetscOptionsGetBool(NULL,prefix,"-pc_asm_print_subdomains",&flg,NULL);CHKERRQ(ierr); 26192bb6962SLisandro Dalcin if (flg) { ierr = PCASMPrintSubdomains(pc);CHKERRQ(ierr); } 2624b9ad928SBarry Smith 2633d03bb48SJed Brown if (osm->overlap > 0) { 2644b9ad928SBarry Smith /* Extend the "overlapping" regions by a number of steps */ 26592bb6962SLisandro Dalcin ierr = MatIncreaseOverlap(pc->pmat,osm->n_local_true,osm->is,osm->overlap);CHKERRQ(ierr); 2663d03bb48SJed Brown } 2676ed231c7SMatthew Knepley if (osm->sort_indices) { 26892bb6962SLisandro Dalcin for (i=0; i<osm->n_local_true; i++) { 2694b9ad928SBarry Smith ierr = ISSort(osm->is[i]);CHKERRQ(ierr); 2702b691e39SMatthew Knepley if (osm->is_local) { 2712b691e39SMatthew Knepley ierr = ISSort(osm->is_local[i]);CHKERRQ(ierr); 2722b691e39SMatthew Knepley } 2734b9ad928SBarry Smith } 2746ed231c7SMatthew Knepley } 2754b9ad928SBarry Smith 276f6991133SBarry Smith if (!osm->ksp) { 27778904715SLisandro Dalcin /* Create the local solvers */ 278785e854fSJed Brown ierr = PetscMalloc1(osm->n_local_true,&osm->ksp);CHKERRQ(ierr); 279feb221f8SDmitry Karpeev if (domain_dm) { 280feb221f8SDmitry Karpeev ierr = PetscInfo(pc,"Setting up ASM subproblems using the embedded DM\n");CHKERRQ(ierr); 281feb221f8SDmitry Karpeev } 28292bb6962SLisandro Dalcin for (i=0; i<osm->n_local_true; i++) { 2834b9ad928SBarry Smith ierr = KSPCreate(PETSC_COMM_SELF,&ksp);CHKERRQ(ierr); 284422a814eSBarry Smith ierr = KSPSetErrorIfNotConverged(ksp,pc->erroriffailure);CHKERRQ(ierr); 2853bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ksp);CHKERRQ(ierr); 28692bb6962SLisandro Dalcin ierr = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,1);CHKERRQ(ierr); 2874b9ad928SBarry Smith ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr); 2884b9ad928SBarry Smith ierr = KSPGetPC(ksp,&subpc);CHKERRQ(ierr); 2894b9ad928SBarry Smith ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 2904b9ad928SBarry Smith ierr = KSPSetOptionsPrefix(ksp,prefix);CHKERRQ(ierr); 2914b9ad928SBarry Smith ierr = KSPAppendOptionsPrefix(ksp,"sub_");CHKERRQ(ierr); 292feb221f8SDmitry Karpeev if (domain_dm) { 293feb221f8SDmitry Karpeev ierr = KSPSetDM(ksp, domain_dm[i]);CHKERRQ(ierr); 294feb221f8SDmitry Karpeev ierr = KSPSetDMActive(ksp, PETSC_FALSE);CHKERRQ(ierr); 295feb221f8SDmitry Karpeev ierr = DMDestroy(&domain_dm[i]);CHKERRQ(ierr); 296feb221f8SDmitry Karpeev } 2974b9ad928SBarry Smith osm->ksp[i] = ksp; 2984b9ad928SBarry Smith } 299feb221f8SDmitry Karpeev if (domain_dm) { 300feb221f8SDmitry Karpeev ierr = PetscFree(domain_dm);CHKERRQ(ierr); 301feb221f8SDmitry Karpeev } 302f6991133SBarry Smith } 3031dd8081eSeaulisa 304fb745f2cSMatthew G. Knepley ierr = ISConcatenate(PETSC_COMM_SELF, osm->n_local_true, osm->is, &osm->lis);CHKERRQ(ierr); 305fb745f2cSMatthew G. Knepley ierr = ISSortRemoveDups(osm->lis);CHKERRQ(ierr); 306fb745f2cSMatthew G. Knepley ierr = ISGetLocalSize(osm->lis, &m);CHKERRQ(ierr); 307fb745f2cSMatthew G. Knepley ierr = VecCreateSeq(PETSC_COMM_SELF, m, &osm->lx);CHKERRQ(ierr); 308fb745f2cSMatthew G. Knepley ierr = VecDuplicate(osm->lx, &osm->ly);CHKERRQ(ierr); 3091dd8081eSeaulisa 3104b9ad928SBarry Smith scall = MAT_INITIAL_MATRIX; 3114b9ad928SBarry Smith } else { 3124b9ad928SBarry Smith /* 3134b9ad928SBarry Smith Destroy the blocks from the previous iteration 3144b9ad928SBarry Smith */ 315e09e08ccSBarry Smith if (pc->flag == DIFFERENT_NONZERO_PATTERN) { 3164b9ad928SBarry Smith ierr = MatDestroyMatrices(osm->n_local_true,&osm->pmat);CHKERRQ(ierr); 3174b9ad928SBarry Smith scall = MAT_INITIAL_MATRIX; 3184b9ad928SBarry Smith } 3194b9ad928SBarry Smith } 3204b9ad928SBarry Smith 32192bb6962SLisandro Dalcin /* 32292bb6962SLisandro Dalcin Extract out the submatrices 32392bb6962SLisandro Dalcin */ 3247dae84e0SHong Zhang ierr = MatCreateSubMatrices(pc->pmat,osm->n_local_true,osm->is,osm->is,scall,&osm->pmat);CHKERRQ(ierr); 32592bb6962SLisandro Dalcin if (scall == MAT_INITIAL_MATRIX) { 32692bb6962SLisandro Dalcin ierr = PetscObjectGetOptionsPrefix((PetscObject)pc->pmat,&pprefix);CHKERRQ(ierr); 32792bb6962SLisandro Dalcin for (i=0; i<osm->n_local_true; i++) { 3283bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)osm->pmat[i]);CHKERRQ(ierr); 32978904715SLisandro Dalcin ierr = PetscObjectSetOptionsPrefix((PetscObject)osm->pmat[i],pprefix);CHKERRQ(ierr); 33092bb6962SLisandro Dalcin } 33192bb6962SLisandro Dalcin } 33280ec0b7dSPatrick Sanan 33380ec0b7dSPatrick Sanan /* Convert the types of the submatrices (if needbe) */ 33480ec0b7dSPatrick Sanan if (osm->sub_mat_type) { 33580ec0b7dSPatrick Sanan for (i=0; i<osm->n_local_true; i++) { 33680ec0b7dSPatrick Sanan ierr = MatConvert(osm->pmat[i],osm->sub_mat_type,MAT_INPLACE_MATRIX,&(osm->pmat[i]));CHKERRQ(ierr); 33780ec0b7dSPatrick Sanan } 33880ec0b7dSPatrick Sanan } 33980ec0b7dSPatrick Sanan 34080ec0b7dSPatrick Sanan if(!pc->setupcalled){ 34180ec0b7dSPatrick Sanan /* Create the local work vectors (from the local matrices) and scatter contexts */ 34280ec0b7dSPatrick Sanan ierr = MatCreateVecs(pc->pmat,&vec,0);CHKERRQ(ierr); 3435b3c0d42Seaulisa 3441dd8081eSeaulisa 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()"); 3451dd8081eSeaulisa if (osm->is_local && osm->type == PC_ASM_RESTRICT && osm->loctype == PC_COMPOSITE_ADDITIVE) { 3461dd8081eSeaulisa ierr = PetscMalloc1(osm->n_local_true,&osm->lprolongation);CHKERRQ(ierr); 3471dd8081eSeaulisa } 3481dd8081eSeaulisa ierr = PetscMalloc1(osm->n_local_true,&osm->lrestriction);CHKERRQ(ierr); 3491dd8081eSeaulisa ierr = PetscMalloc1(osm->n_local_true,&osm->x);CHKERRQ(ierr); 3501dd8081eSeaulisa ierr = PetscMalloc1(osm->n_local_true,&osm->y);CHKERRQ(ierr); 351347574c9Seaulisa 352347574c9Seaulisa ierr = ISGetLocalSize(osm->lis,&m);CHKERRQ(ierr); 353347574c9Seaulisa ierr = ISCreateStride(PETSC_COMM_SELF,m,0,1,&isl);CHKERRQ(ierr); 35435928de7SBarry Smith ierr = VecScatterCreateWithData(vec,osm->lis,osm->lx,isl,&osm->restriction);CHKERRQ(ierr); 355347574c9Seaulisa ierr = ISDestroy(&isl);CHKERRQ(ierr); 356347574c9Seaulisa 357347574c9Seaulisa 35880ec0b7dSPatrick Sanan for (i=0; i<osm->n_local_true; ++i) { 3595b3c0d42Seaulisa ISLocalToGlobalMapping ltog; 3605b3c0d42Seaulisa IS isll; 3615b3c0d42Seaulisa const PetscInt *idx_is; 3625b3c0d42Seaulisa PetscInt *idx_lis,nout; 3635b3c0d42Seaulisa 36455817e79SBarry Smith ierr = ISGetLocalSize(osm->is[i],&m);CHKERRQ(ierr); 36555817e79SBarry Smith ierr = MatCreateVecs(osm->pmat[i],&osm->x[i],NULL);CHKERRQ(ierr); 36655817e79SBarry Smith ierr = VecDuplicate(osm->x[i],&osm->y[i]);CHKERRQ(ierr); 36755817e79SBarry Smith 368b0de9f23SBarry Smith /* generate a scatter from ly to y[i] picking all the overlapping is[i] entries */ 3695b3c0d42Seaulisa ierr = ISLocalToGlobalMappingCreateIS(osm->lis,<og);CHKERRQ(ierr); 3705b3c0d42Seaulisa ierr = ISGetLocalSize(osm->is[i],&m);CHKERRQ(ierr); 3715b3c0d42Seaulisa ierr = ISGetIndices(osm->is[i], &idx_is);CHKERRQ(ierr); 3725b3c0d42Seaulisa ierr = PetscMalloc1(m,&idx_lis);CHKERRQ(ierr); 3735b3c0d42Seaulisa ierr = ISGlobalToLocalMappingApply(ltog,IS_GTOLM_DROP,m,idx_is,&nout,idx_lis);CHKERRQ(ierr); 3745b3c0d42Seaulisa if (nout != m) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"is not a subset of lis"); 3755b3c0d42Seaulisa ierr = ISRestoreIndices(osm->is[i], &idx_is);CHKERRQ(ierr); 3765b3c0d42Seaulisa ierr = ISCreateGeneral(PETSC_COMM_SELF,m,idx_lis,PETSC_OWN_POINTER,&isll);CHKERRQ(ierr); 377d8b3b5e3Seaulisa ierr = ISLocalToGlobalMappingDestroy(<og);CHKERRQ(ierr); 3785b3c0d42Seaulisa ierr = ISCreateStride(PETSC_COMM_SELF,m,0,1,&isl);CHKERRQ(ierr); 37935928de7SBarry Smith ierr = VecScatterCreateWithData(osm->ly,isll,osm->y[i],isl,&osm->lrestriction[i]);CHKERRQ(ierr); 3805b3c0d42Seaulisa ierr = ISDestroy(&isll);CHKERRQ(ierr); 3815b3c0d42Seaulisa ierr = ISDestroy(&isl);CHKERRQ(ierr); 382b0de9f23SBarry Smith if (osm->lprolongation) { /* generate a scatter from y[i] to ly picking only the the non-overalapping is_local[i] entries */ 383d8b3b5e3Seaulisa ISLocalToGlobalMapping ltog; 384d8b3b5e3Seaulisa IS isll,isll_local; 385d8b3b5e3Seaulisa const PetscInt *idx_local; 386d8b3b5e3Seaulisa PetscInt *idx1, *idx2, nout; 387d8b3b5e3Seaulisa 388d8b3b5e3Seaulisa ierr = ISGetLocalSize(osm->is_local[i],&m_local);CHKERRQ(ierr); 389d8b3b5e3Seaulisa ierr = ISGetIndices(osm->is_local[i], &idx_local);CHKERRQ(ierr); 390d8b3b5e3Seaulisa 391d8b3b5e3Seaulisa ierr = ISLocalToGlobalMappingCreateIS(osm->is[i],<og);CHKERRQ(ierr); 392d8b3b5e3Seaulisa ierr = PetscMalloc1(m_local,&idx1);CHKERRQ(ierr); 393d8b3b5e3Seaulisa ierr = ISGlobalToLocalMappingApply(ltog,IS_GTOLM_DROP,m_local,idx_local,&nout,idx1);CHKERRQ(ierr); 394d8b3b5e3Seaulisa ierr = ISLocalToGlobalMappingDestroy(<og);CHKERRQ(ierr); 395d8b3b5e3Seaulisa if (nout != m_local) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"is_local not a subset of is"); 396d8b3b5e3Seaulisa ierr = ISCreateGeneral(PETSC_COMM_SELF,m_local,idx1,PETSC_OWN_POINTER,&isll);CHKERRQ(ierr); 397d8b3b5e3Seaulisa 398d8b3b5e3Seaulisa ierr = ISLocalToGlobalMappingCreateIS(osm->lis,<og);CHKERRQ(ierr); 399d8b3b5e3Seaulisa ierr = PetscMalloc1(m_local,&idx2);CHKERRQ(ierr); 400d8b3b5e3Seaulisa ierr = ISGlobalToLocalMappingApply(ltog,IS_GTOLM_DROP,m_local,idx_local,&nout,idx2);CHKERRQ(ierr); 401d8b3b5e3Seaulisa ierr = ISLocalToGlobalMappingDestroy(<og);CHKERRQ(ierr); 402d8b3b5e3Seaulisa if (nout != m_local) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"is_local not a subset of lis"); 403d8b3b5e3Seaulisa ierr = ISCreateGeneral(PETSC_COMM_SELF,m_local,idx2,PETSC_OWN_POINTER,&isll_local);CHKERRQ(ierr); 404d8b3b5e3Seaulisa 405d8b3b5e3Seaulisa ierr = ISRestoreIndices(osm->is_local[i], &idx_local);CHKERRQ(ierr); 40635928de7SBarry Smith ierr = VecScatterCreateWithData(osm->y[i],isll,osm->ly,isll_local,&osm->lprolongation[i]);CHKERRQ(ierr); 407d8b3b5e3Seaulisa 408d8b3b5e3Seaulisa ierr = ISDestroy(&isll);CHKERRQ(ierr); 409d8b3b5e3Seaulisa ierr = ISDestroy(&isll_local);CHKERRQ(ierr); 410d8b3b5e3Seaulisa } 41180ec0b7dSPatrick Sanan } 41280ec0b7dSPatrick Sanan ierr = VecDestroy(&vec);CHKERRQ(ierr); 41380ec0b7dSPatrick Sanan } 41480ec0b7dSPatrick Sanan 415fb745f2cSMatthew G. Knepley if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE) { 416235cc4ceSMatthew G. Knepley IS *cis; 417235cc4ceSMatthew G. Knepley PetscInt c; 418235cc4ceSMatthew G. Knepley 419235cc4ceSMatthew G. Knepley ierr = PetscMalloc1(osm->n_local_true, &cis);CHKERRQ(ierr); 420235cc4ceSMatthew G. Knepley for (c = 0; c < osm->n_local_true; ++c) cis[c] = osm->lis; 4217dae84e0SHong Zhang ierr = MatCreateSubMatrices(pc->pmat, osm->n_local_true, osm->is, cis, scall, &osm->lmats);CHKERRQ(ierr); 422235cc4ceSMatthew G. Knepley ierr = PetscFree(cis);CHKERRQ(ierr); 423fb745f2cSMatthew G. Knepley } 4244b9ad928SBarry Smith 4254b9ad928SBarry Smith /* Return control to the user so that the submatrices can be modified (e.g., to apply 4264b9ad928SBarry Smith different boundary conditions for the submatrices than for the global problem) */ 42792bb6962SLisandro Dalcin ierr = PCModifySubMatrices(pc,osm->n_local_true,osm->is,osm->is,osm->pmat,pc->modifysubmatricesP);CHKERRQ(ierr); 4284b9ad928SBarry Smith 42992bb6962SLisandro Dalcin /* 43092bb6962SLisandro Dalcin Loop over subdomains putting them into local ksp 43192bb6962SLisandro Dalcin */ 43292bb6962SLisandro Dalcin for (i=0; i<osm->n_local_true; i++) { 43323ee1639SBarry Smith ierr = KSPSetOperators(osm->ksp[i],osm->pmat[i],osm->pmat[i]);CHKERRQ(ierr); 43492bb6962SLisandro Dalcin if (!pc->setupcalled) { 435bf108f30SBarry Smith ierr = KSPSetFromOptions(osm->ksp[i]);CHKERRQ(ierr); 4364b9ad928SBarry Smith } 43792bb6962SLisandro Dalcin } 4384b9ad928SBarry Smith PetscFunctionReturn(0); 4394b9ad928SBarry Smith } 4404b9ad928SBarry Smith 4416849ba73SBarry Smith static PetscErrorCode PCSetUpOnBlocks_ASM(PC pc) 4424b9ad928SBarry Smith { 4434b9ad928SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 4446849ba73SBarry Smith PetscErrorCode ierr; 44513f74950SBarry Smith PetscInt i; 446539c167fSBarry Smith KSPConvergedReason reason; 4474b9ad928SBarry Smith 4484b9ad928SBarry Smith PetscFunctionBegin; 4494b9ad928SBarry Smith for (i=0; i<osm->n_local_true; i++) { 4504b9ad928SBarry Smith ierr = KSPSetUp(osm->ksp[i]);CHKERRQ(ierr); 451539c167fSBarry Smith ierr = KSPGetConvergedReason(osm->ksp[i],&reason);CHKERRQ(ierr); 452*c0decd05SBarry Smith if (reason == KSP_DIVERGED_PC_FAILED) { 453261222b2SHong Zhang pc->failedreason = PC_SUBPC_ERROR; 454e0eafd54SHong Zhang } 4554b9ad928SBarry Smith } 4564b9ad928SBarry Smith PetscFunctionReturn(0); 4574b9ad928SBarry Smith } 4584b9ad928SBarry Smith 4596849ba73SBarry Smith static PetscErrorCode PCApply_ASM(PC pc,Vec x,Vec y) 4604b9ad928SBarry Smith { 4614b9ad928SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 4626849ba73SBarry Smith PetscErrorCode ierr; 4631dd8081eSeaulisa PetscInt i,n_local_true = osm->n_local_true; 4644b9ad928SBarry Smith ScatterMode forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE; 4654b9ad928SBarry Smith 4664b9ad928SBarry Smith PetscFunctionBegin; 4674b9ad928SBarry Smith /* 4684b9ad928SBarry Smith Support for limiting the restriction or interpolation to only local 4694b9ad928SBarry Smith subdomain values (leaving the other values 0). 4704b9ad928SBarry Smith */ 4714b9ad928SBarry Smith if (!(osm->type & PC_ASM_RESTRICT)) { 4724b9ad928SBarry Smith forward = SCATTER_FORWARD_LOCAL; 4734b9ad928SBarry Smith /* have to zero the work RHS since scatter may leave some slots empty */ 4741dd8081eSeaulisa ierr = VecSet(osm->lx, 0.0);CHKERRQ(ierr); 4754b9ad928SBarry Smith } 476347574c9Seaulisa if (!(osm->type & PC_ASM_INTERPOLATE)) { 477347574c9Seaulisa reverse = SCATTER_REVERSE_LOCAL; 478347574c9Seaulisa } 4794b9ad928SBarry Smith 480347574c9Seaulisa if(osm->loctype == PC_COMPOSITE_MULTIPLICATIVE || osm->loctype == PC_COMPOSITE_ADDITIVE){ 481b0de9f23SBarry Smith /* zero the global and the local solutions */ 48212cd4985SMatthew G. Knepley ierr = VecZeroEntries(y);CHKERRQ(ierr); 4835b3c0d42Seaulisa ierr = VecSet(osm->ly, 0.0);CHKERRQ(ierr); 484347574c9Seaulisa 485b0de9f23SBarry Smith /* Copy the global RHS to local RHS including the ghost nodes */ 4861dd8081eSeaulisa ierr = VecScatterBegin(osm->restriction, x, osm->lx, INSERT_VALUES, forward);CHKERRQ(ierr); 4871dd8081eSeaulisa ierr = VecScatterEnd(osm->restriction, x, osm->lx, INSERT_VALUES, forward);CHKERRQ(ierr); 488347574c9Seaulisa 489b0de9f23SBarry Smith /* Restrict local RHS to the overlapping 0-block RHS */ 4901dd8081eSeaulisa ierr = VecScatterBegin(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward);CHKERRQ(ierr); 4911dd8081eSeaulisa ierr = VecScatterEnd(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward);CHKERRQ(ierr); 492d8b3b5e3Seaulisa 49312cd4985SMatthew G. Knepley /* do the local solves */ 49412cd4985SMatthew G. Knepley for (i = 0; i < n_local_true; ++i) { 495347574c9Seaulisa 496b0de9f23SBarry Smith /* solve the overlapping i-block */ 497fa2bb9feSLisandro Dalcin ierr = PetscLogEventBegin(PC_ApplyOnBlocks,osm->ksp[i],osm->x[i],osm->y[i],0);CHKERRQ(ierr); 49812cd4985SMatthew G. Knepley ierr = KSPSolve(osm->ksp[i], osm->x[i], osm->y[i]);CHKERRQ(ierr); 499*c0decd05SBarry Smith ierr = KSPCheckSolve(osm->ksp[i],pc,osm->y[i]);CHKERRQ(ierr); 500fa2bb9feSLisandro Dalcin ierr = PetscLogEventEnd(PC_ApplyOnBlocks,osm->ksp[i],osm->x[i],osm->y[i],0);CHKERRQ(ierr); 501d8b3b5e3Seaulisa 502b0de9f23SBarry Smith if (osm->lprolongation) { /* interpolate the non-overalapping i-block solution to the local solution (only for restrictive additive) */ 5031dd8081eSeaulisa ierr = VecScatterBegin(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward);CHKERRQ(ierr); 5041dd8081eSeaulisa ierr = VecScatterEnd(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward);CHKERRQ(ierr); 505d8b3b5e3Seaulisa } 506b0de9f23SBarry Smith else{ /* interpolate the overalapping i-block solution to the local solution */ 5071dd8081eSeaulisa ierr = VecScatterBegin(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse);CHKERRQ(ierr); 5081dd8081eSeaulisa ierr = VecScatterEnd(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse);CHKERRQ(ierr); 509d8b3b5e3Seaulisa } 510347574c9Seaulisa 511347574c9Seaulisa if (i < n_local_true-1) { 512b0de9f23SBarry Smith /* Restrict local RHS to the overlapping (i+1)-block RHS */ 5131dd8081eSeaulisa ierr = VecScatterBegin(osm->lrestriction[i+1], osm->lx, osm->x[i+1], INSERT_VALUES, forward);CHKERRQ(ierr); 5141dd8081eSeaulisa ierr = VecScatterEnd(osm->lrestriction[i+1], osm->lx, osm->x[i+1], INSERT_VALUES, forward);CHKERRQ(ierr); 515347574c9Seaulisa 516347574c9Seaulisa if ( osm->loctype == PC_COMPOSITE_MULTIPLICATIVE){ 517b0de9f23SBarry Smith /* udpdate the overlapping (i+1)-block RHS using the current local solution */ 518347574c9Seaulisa ierr = MatMult(osm->lmats[i+1], osm->ly, osm->y[i+1]);CHKERRQ(ierr); 519347574c9Seaulisa ierr = VecAXPBY(osm->x[i+1],-1.,1., osm->y[i+1]); CHKERRQ(ierr); 5207c3d802fSMatthew G. Knepley } 52112cd4985SMatthew G. Knepley } 52212cd4985SMatthew G. Knepley } 523b0de9f23SBarry Smith /* Add the local solution to the global solution including the ghost nodes */ 5241dd8081eSeaulisa ierr = VecScatterBegin(osm->restriction, osm->ly, y, ADD_VALUES, reverse);CHKERRQ(ierr); 5251dd8081eSeaulisa ierr = VecScatterEnd(osm->restriction, osm->ly, y, ADD_VALUES, reverse);CHKERRQ(ierr); 526347574c9Seaulisa }else{ 527347574c9Seaulisa SETERRQ1(PetscObjectComm((PetscObject) pc), PETSC_ERR_ARG_WRONG, "Invalid local composition type: %s", PCCompositeTypes[osm->loctype]); 52812cd4985SMatthew G. Knepley } 5294b9ad928SBarry Smith PetscFunctionReturn(0); 5304b9ad928SBarry Smith } 5314b9ad928SBarry Smith 5326849ba73SBarry Smith static PetscErrorCode PCApplyTranspose_ASM(PC pc,Vec x,Vec y) 5334b9ad928SBarry Smith { 5344b9ad928SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 5356849ba73SBarry Smith PetscErrorCode ierr; 5361dd8081eSeaulisa PetscInt i,n_local_true = osm->n_local_true; 5374b9ad928SBarry Smith ScatterMode forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE; 5384b9ad928SBarry Smith 5394b9ad928SBarry Smith PetscFunctionBegin; 5404b9ad928SBarry Smith /* 5414b9ad928SBarry Smith Support for limiting the restriction or interpolation to only local 5424b9ad928SBarry Smith subdomain values (leaving the other values 0). 5434b9ad928SBarry Smith 5444b9ad928SBarry Smith Note: these are reversed from the PCApply_ASM() because we are applying the 5454b9ad928SBarry Smith transpose of the three terms 5464b9ad928SBarry Smith */ 547d8b3b5e3Seaulisa 5484b9ad928SBarry Smith if (!(osm->type & PC_ASM_INTERPOLATE)) { 5494b9ad928SBarry Smith forward = SCATTER_FORWARD_LOCAL; 5504b9ad928SBarry Smith /* have to zero the work RHS since scatter may leave some slots empty */ 5511dd8081eSeaulisa ierr = VecSet(osm->lx, 0.0);CHKERRQ(ierr); 5524b9ad928SBarry Smith } 5532fa5cd67SKarl Rupp if (!(osm->type & PC_ASM_RESTRICT)) reverse = SCATTER_REVERSE_LOCAL; 5544b9ad928SBarry Smith 555b0de9f23SBarry Smith /* zero the global and the local solutions */ 55610bd88b9SJed Brown ierr = VecZeroEntries(y);CHKERRQ(ierr); 557d8b3b5e3Seaulisa ierr = VecSet(osm->ly, 0.0);CHKERRQ(ierr); 558d8b3b5e3Seaulisa 559b0de9f23SBarry Smith /* Copy the global RHS to local RHS including the ghost nodes */ 5601dd8081eSeaulisa ierr = VecScatterBegin(osm->restriction, x, osm->lx, INSERT_VALUES, forward);CHKERRQ(ierr); 5611dd8081eSeaulisa ierr = VecScatterEnd(osm->restriction, x, osm->lx, INSERT_VALUES, forward);CHKERRQ(ierr); 562d8b3b5e3Seaulisa 563b0de9f23SBarry Smith /* Restrict local RHS to the overlapping 0-block RHS */ 5641dd8081eSeaulisa ierr = VecScatterBegin(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward);CHKERRQ(ierr); 5651dd8081eSeaulisa ierr = VecScatterEnd(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward);CHKERRQ(ierr); 566d8b3b5e3Seaulisa 5674b9ad928SBarry Smith /* do the local solves */ 568d8b3b5e3Seaulisa for (i = 0; i < n_local_true; ++i) { 569d8b3b5e3Seaulisa 570b0de9f23SBarry Smith /* solve the overlapping i-block */ 571fa2bb9feSLisandro Dalcin ierr = PetscLogEventBegin(PC_ApplyOnBlocks,osm->ksp[i],osm->x[i],osm->y[i],0);CHKERRQ(ierr); 572e1bcd54cSBarry Smith ierr = KSPSolveTranspose(osm->ksp[i], osm->x[i], osm->y[i]);CHKERRQ(ierr); 573*c0decd05SBarry Smith ierr = KSPCheckSolve(osm->ksp[i],pc,osm->y[i]);CHKERRQ(ierr); 574fa2bb9feSLisandro Dalcin ierr = PetscLogEventEnd(PC_ApplyOnBlocks,osm->ksp[i],osm->x[i],osm->y[i],0);CHKERRQ(ierr); 575d8b3b5e3Seaulisa 576b0de9f23SBarry Smith if (osm->lprolongation) { /* interpolate the non-overalapping i-block solution to the local solution */ 5771dd8081eSeaulisa ierr = VecScatterBegin(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward);CHKERRQ(ierr); 5781dd8081eSeaulisa ierr = VecScatterEnd(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward);CHKERRQ(ierr); 579b9845e0eSMatthew Knepley } 580b0de9f23SBarry Smith else{ /* interpolate the overalapping i-block solution to the local solution */ 5811dd8081eSeaulisa ierr = VecScatterBegin(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse);CHKERRQ(ierr); 5821dd8081eSeaulisa ierr = VecScatterEnd(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse);CHKERRQ(ierr); 5834b9ad928SBarry Smith } 584d8b3b5e3Seaulisa 585d8b3b5e3Seaulisa if (i < n_local_true-1) { 586b0de9f23SBarry Smith /* Restrict local RHS to the overlapping (i+1)-block RHS */ 5871dd8081eSeaulisa ierr = VecScatterBegin(osm->lrestriction[i+1], osm->lx, osm->x[i+1], INSERT_VALUES, forward);CHKERRQ(ierr); 5881dd8081eSeaulisa ierr = VecScatterEnd(osm->lrestriction[i+1], osm->lx, osm->x[i+1], INSERT_VALUES, forward);CHKERRQ(ierr); 5894b9ad928SBarry Smith } 5904b9ad928SBarry Smith } 591b0de9f23SBarry Smith /* Add the local solution to the global solution including the ghost nodes */ 5921dd8081eSeaulisa ierr = VecScatterBegin(osm->restriction, osm->ly, y, ADD_VALUES, reverse);CHKERRQ(ierr); 5931dd8081eSeaulisa ierr = VecScatterEnd(osm->restriction, osm->ly, y, ADD_VALUES, reverse);CHKERRQ(ierr); 594d8b3b5e3Seaulisa 5954b9ad928SBarry Smith PetscFunctionReturn(0); 596d8b3b5e3Seaulisa 5974b9ad928SBarry Smith } 5984b9ad928SBarry Smith 599e91c6855SBarry Smith static PetscErrorCode PCReset_ASM(PC pc) 6004b9ad928SBarry Smith { 6014b9ad928SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 6026849ba73SBarry Smith PetscErrorCode ierr; 60313f74950SBarry Smith PetscInt i; 6044b9ad928SBarry Smith 6054b9ad928SBarry Smith PetscFunctionBegin; 60692bb6962SLisandro Dalcin if (osm->ksp) { 60792bb6962SLisandro Dalcin for (i=0; i<osm->n_local_true; i++) { 608e91c6855SBarry Smith ierr = KSPReset(osm->ksp[i]);CHKERRQ(ierr); 60992bb6962SLisandro Dalcin } 61092bb6962SLisandro Dalcin } 611e09e08ccSBarry Smith if (osm->pmat) { 61292bb6962SLisandro Dalcin if (osm->n_local_true > 0) { 61330a70a9aSHong Zhang ierr = MatDestroySubMatrices(osm->n_local_true,&osm->pmat);CHKERRQ(ierr); 61492bb6962SLisandro Dalcin } 61592bb6962SLisandro Dalcin } 6161dd8081eSeaulisa if (osm->lrestriction) { 6171dd8081eSeaulisa ierr = VecScatterDestroy(&osm->restriction);CHKERRQ(ierr); 6181dd8081eSeaulisa for (i=0; i<osm->n_local_true; i++) { 6191dd8081eSeaulisa ierr = VecScatterDestroy(&osm->lrestriction[i]);CHKERRQ(ierr); 6201dd8081eSeaulisa if (osm->lprolongation) {ierr = VecScatterDestroy(&osm->lprolongation[i]);CHKERRQ(ierr);} 621fcfd50ebSBarry Smith ierr = VecDestroy(&osm->x[i]);CHKERRQ(ierr); 622fcfd50ebSBarry Smith ierr = VecDestroy(&osm->y[i]);CHKERRQ(ierr); 6234b9ad928SBarry Smith } 6241dd8081eSeaulisa ierr = PetscFree(osm->lrestriction);CHKERRQ(ierr); 6251dd8081eSeaulisa if (osm->lprolongation) {ierr = PetscFree(osm->lprolongation);CHKERRQ(ierr);} 62605b42c5fSBarry Smith ierr = PetscFree(osm->x);CHKERRQ(ierr); 62778904715SLisandro Dalcin ierr = PetscFree(osm->y);CHKERRQ(ierr); 6281dd8081eSeaulisa 62992bb6962SLisandro Dalcin } 6302b691e39SMatthew Knepley ierr = PCASMDestroySubdomains(osm->n_local_true,osm->is,osm->is_local);CHKERRQ(ierr); 631fb745f2cSMatthew G. Knepley ierr = ISDestroy(&osm->lis);CHKERRQ(ierr); 632fb745f2cSMatthew G. Knepley ierr = VecDestroy(&osm->lx);CHKERRQ(ierr); 633fb745f2cSMatthew G. Knepley ierr = VecDestroy(&osm->ly);CHKERRQ(ierr); 634347574c9Seaulisa if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE) { 635347574c9Seaulisa ierr = MatDestroyMatrices(osm->n_local_true, &osm->lmats);CHKERRQ(ierr); 636fb745f2cSMatthew G. Knepley } 6372fa5cd67SKarl Rupp 63880ec0b7dSPatrick Sanan ierr = PetscFree(osm->sub_mat_type);CHKERRQ(ierr); 63980ec0b7dSPatrick Sanan 640e91c6855SBarry Smith osm->is = 0; 641e91c6855SBarry Smith osm->is_local = 0; 642e91c6855SBarry Smith PetscFunctionReturn(0); 643e91c6855SBarry Smith } 644e91c6855SBarry Smith 645e91c6855SBarry Smith static PetscErrorCode PCDestroy_ASM(PC pc) 646e91c6855SBarry Smith { 647e91c6855SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 648e91c6855SBarry Smith PetscErrorCode ierr; 649e91c6855SBarry Smith PetscInt i; 650e91c6855SBarry Smith 651e91c6855SBarry Smith PetscFunctionBegin; 652e91c6855SBarry Smith ierr = PCReset_ASM(pc);CHKERRQ(ierr); 653e91c6855SBarry Smith if (osm->ksp) { 654e91c6855SBarry Smith for (i=0; i<osm->n_local_true; i++) { 655fcfd50ebSBarry Smith ierr = KSPDestroy(&osm->ksp[i]);CHKERRQ(ierr); 656e91c6855SBarry Smith } 657e91c6855SBarry Smith ierr = PetscFree(osm->ksp);CHKERRQ(ierr); 658e91c6855SBarry Smith } 659e91c6855SBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 6604b9ad928SBarry Smith PetscFunctionReturn(0); 6614b9ad928SBarry Smith } 6624b9ad928SBarry Smith 6634416b707SBarry Smith static PetscErrorCode PCSetFromOptions_ASM(PetscOptionItems *PetscOptionsObject,PC pc) 6644b9ad928SBarry Smith { 6654b9ad928SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 6666849ba73SBarry Smith PetscErrorCode ierr; 6679dcbbd2bSBarry Smith PetscInt blocks,ovl; 668ace3abfcSBarry Smith PetscBool symset,flg; 66992bb6962SLisandro Dalcin PCASMType asmtype; 67012cd4985SMatthew G. Knepley PCCompositeType loctype; 67180ec0b7dSPatrick Sanan char sub_mat_type[256]; 6724b9ad928SBarry Smith 6734b9ad928SBarry Smith PetscFunctionBegin; 674bf108f30SBarry Smith /* set the type to symmetric if matrix is symmetric */ 67592bb6962SLisandro Dalcin if (!osm->type_set && pc->pmat) { 67692bb6962SLisandro Dalcin ierr = MatIsSymmetricKnown(pc->pmat,&symset,&flg);CHKERRQ(ierr); 6772fa5cd67SKarl Rupp if (symset && flg) osm->type = PC_ASM_BASIC; 678bf108f30SBarry Smith } 679e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"Additive Schwarz options");CHKERRQ(ierr); 680d709ea83SDmitry Karpeev ierr = PetscOptionsBool("-pc_asm_dm_subdomains","Use DMCreateDomainDecomposition() to define subdomains","PCASMSetDMSubdomains",osm->dm_subdomains,&osm->dm_subdomains,&flg);CHKERRQ(ierr); 68190d69ab7SBarry Smith ierr = PetscOptionsInt("-pc_asm_blocks","Number of subdomains","PCASMSetTotalSubdomains",osm->n,&blocks,&flg);CHKERRQ(ierr); 68265db9045SDmitry Karpeev if (flg) { 683f77a5242SKarl Rupp ierr = PCASMSetTotalSubdomains(pc,blocks,NULL,NULL);CHKERRQ(ierr); 684d709ea83SDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 68565db9045SDmitry Karpeev } 686342c94f9SMatthew G. Knepley ierr = PetscOptionsInt("-pc_asm_local_blocks","Number of local subdomains","PCASMSetLocalSubdomains",osm->n_local_true,&blocks,&flg);CHKERRQ(ierr); 687342c94f9SMatthew G. Knepley if (flg) { 688342c94f9SMatthew G. Knepley ierr = PCASMSetLocalSubdomains(pc,blocks,NULL,NULL);CHKERRQ(ierr); 689342c94f9SMatthew G. Knepley osm->dm_subdomains = PETSC_FALSE; 690342c94f9SMatthew G. Knepley } 69190d69ab7SBarry Smith ierr = PetscOptionsInt("-pc_asm_overlap","Number of grid points overlap","PCASMSetOverlap",osm->overlap,&ovl,&flg);CHKERRQ(ierr); 69265db9045SDmitry Karpeev if (flg) { 69365db9045SDmitry Karpeev ierr = PCASMSetOverlap(pc,ovl);CHKERRQ(ierr); 694d709ea83SDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 69565db9045SDmitry Karpeev } 69690d69ab7SBarry Smith flg = PETSC_FALSE; 69790d69ab7SBarry Smith ierr = PetscOptionsEnum("-pc_asm_type","Type of restriction/extension","PCASMSetType",PCASMTypes,(PetscEnum)osm->type,(PetscEnum*)&asmtype,&flg);CHKERRQ(ierr); 69892bb6962SLisandro Dalcin if (flg) {ierr = PCASMSetType(pc,asmtype);CHKERRQ(ierr); } 69912cd4985SMatthew G. Knepley flg = PETSC_FALSE; 70012cd4985SMatthew G. Knepley ierr = PetscOptionsEnum("-pc_asm_local_type","Type of local solver composition","PCASMSetLocalType",PCCompositeTypes,(PetscEnum)osm->loctype,(PetscEnum*)&loctype,&flg);CHKERRQ(ierr); 70112cd4985SMatthew G. Knepley if (flg) {ierr = PCASMSetLocalType(pc,loctype);CHKERRQ(ierr); } 70280ec0b7dSPatrick Sanan ierr = PetscOptionsFList("-pc_asm_sub_mat_type","Subsolve Matrix Type","PCASMSetSubMatType",MatList,NULL,sub_mat_type,256,&flg);CHKERRQ(ierr); 70380ec0b7dSPatrick Sanan if(flg){ 704459726d8SSatish Balay ierr = PCASMSetSubMatType(pc,sub_mat_type);CHKERRQ(ierr); 70580ec0b7dSPatrick Sanan } 7064b9ad928SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 7074b9ad928SBarry Smith PetscFunctionReturn(0); 7084b9ad928SBarry Smith } 7094b9ad928SBarry Smith 7104b9ad928SBarry Smith /*------------------------------------------------------------------------------------*/ 7114b9ad928SBarry Smith 7121e6b0712SBarry Smith static PetscErrorCode PCASMSetLocalSubdomains_ASM(PC pc,PetscInt n,IS is[],IS is_local[]) 7134b9ad928SBarry Smith { 7144b9ad928SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 71592bb6962SLisandro Dalcin PetscErrorCode ierr; 71692bb6962SLisandro Dalcin PetscInt i; 7174b9ad928SBarry Smith 7184b9ad928SBarry Smith PetscFunctionBegin; 719e32f2f54SBarry Smith if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Each process must have 1 or more blocks, n = %D",n); 720ce94432eSBarry 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()."); 721e7e72b3dSBarry Smith 7224b9ad928SBarry Smith if (!pc->setupcalled) { 72392bb6962SLisandro Dalcin if (is) { 72492bb6962SLisandro Dalcin for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)is[i]);CHKERRQ(ierr);} 72592bb6962SLisandro Dalcin } 726832fc9a5SMatthew Knepley if (is_local) { 727832fc9a5SMatthew Knepley for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)is_local[i]);CHKERRQ(ierr);} 728832fc9a5SMatthew Knepley } 7292b691e39SMatthew Knepley ierr = PCASMDestroySubdomains(osm->n_local_true,osm->is,osm->is_local);CHKERRQ(ierr); 7302fa5cd67SKarl Rupp 7314b9ad928SBarry Smith osm->n_local_true = n; 73292bb6962SLisandro Dalcin osm->is = 0; 7332b691e39SMatthew Knepley osm->is_local = 0; 73492bb6962SLisandro Dalcin if (is) { 735785e854fSJed Brown ierr = PetscMalloc1(n,&osm->is);CHKERRQ(ierr); 7362fa5cd67SKarl Rupp for (i=0; i<n; i++) osm->is[i] = is[i]; 7373d03bb48SJed Brown /* Flag indicating that the user has set overlapping subdomains so PCASM should not increase their size. */ 7383d03bb48SJed Brown osm->overlap = -1; 73992bb6962SLisandro Dalcin } 7402b691e39SMatthew Knepley if (is_local) { 741785e854fSJed Brown ierr = PetscMalloc1(n,&osm->is_local);CHKERRQ(ierr); 7422fa5cd67SKarl Rupp for (i=0; i<n; i++) osm->is_local[i] = is_local[i]; 743a35b7b57SDmitry Karpeev if (!is) { 744785e854fSJed Brown ierr = PetscMalloc1(osm->n_local_true,&osm->is);CHKERRQ(ierr); 745a35b7b57SDmitry Karpeev for (i=0; i<osm->n_local_true; i++) { 746a35b7b57SDmitry Karpeev if (osm->overlap > 0) { /* With positive overlap, osm->is[i] will be modified */ 747a35b7b57SDmitry Karpeev ierr = ISDuplicate(osm->is_local[i],&osm->is[i]);CHKERRQ(ierr); 748a35b7b57SDmitry Karpeev ierr = ISCopy(osm->is_local[i],osm->is[i]);CHKERRQ(ierr); 749a35b7b57SDmitry Karpeev } else { 750a35b7b57SDmitry Karpeev ierr = PetscObjectReference((PetscObject)osm->is_local[i]);CHKERRQ(ierr); 751a35b7b57SDmitry Karpeev osm->is[i] = osm->is_local[i]; 752a35b7b57SDmitry Karpeev } 753a35b7b57SDmitry Karpeev } 754a35b7b57SDmitry Karpeev } 7552b691e39SMatthew Knepley } 7564b9ad928SBarry Smith } 7574b9ad928SBarry Smith PetscFunctionReturn(0); 7584b9ad928SBarry Smith } 7594b9ad928SBarry Smith 7601e6b0712SBarry Smith static PetscErrorCode PCASMSetTotalSubdomains_ASM(PC pc,PetscInt N,IS *is,IS *is_local) 7614b9ad928SBarry Smith { 7624b9ad928SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 7636849ba73SBarry Smith PetscErrorCode ierr; 76413f74950SBarry Smith PetscMPIInt rank,size; 76578904715SLisandro Dalcin PetscInt n; 7664b9ad928SBarry Smith 7674b9ad928SBarry Smith PetscFunctionBegin; 768ce94432eSBarry Smith if (N < 1) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Number of total blocks must be > 0, N = %D",N); 769ce94432eSBarry 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."); 7704b9ad928SBarry Smith 7714b9ad928SBarry Smith /* 772880770e9SJed Brown Split the subdomains equally among all processors 7734b9ad928SBarry Smith */ 774ce94432eSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);CHKERRQ(ierr); 775ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr); 7764b9ad928SBarry Smith n = N/size + ((N % size) > rank); 777e32f2f54SBarry 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); 778e32f2f54SBarry Smith if (pc->setupcalled && n != osm->n_local_true) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"PCASMSetTotalSubdomains() should be called before PCSetUp()."); 7794b9ad928SBarry Smith if (!pc->setupcalled) { 7802b691e39SMatthew Knepley ierr = PCASMDestroySubdomains(osm->n_local_true,osm->is,osm->is_local);CHKERRQ(ierr); 7812fa5cd67SKarl Rupp 7824b9ad928SBarry Smith osm->n_local_true = n; 7834b9ad928SBarry Smith osm->is = 0; 7842b691e39SMatthew Knepley osm->is_local = 0; 7854b9ad928SBarry Smith } 7864b9ad928SBarry Smith PetscFunctionReturn(0); 7874b9ad928SBarry Smith } 7884b9ad928SBarry Smith 7891e6b0712SBarry Smith static PetscErrorCode PCASMSetOverlap_ASM(PC pc,PetscInt ovl) 7904b9ad928SBarry Smith { 79192bb6962SLisandro Dalcin PC_ASM *osm = (PC_ASM*)pc->data; 7924b9ad928SBarry Smith 7934b9ad928SBarry Smith PetscFunctionBegin; 794ce94432eSBarry Smith if (ovl < 0) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap value requested"); 795ce94432eSBarry Smith if (pc->setupcalled && ovl != osm->overlap) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCASMSetOverlap() should be called before PCSetUp()."); 7962fa5cd67SKarl Rupp if (!pc->setupcalled) osm->overlap = ovl; 7974b9ad928SBarry Smith PetscFunctionReturn(0); 7984b9ad928SBarry Smith } 7994b9ad928SBarry Smith 8001e6b0712SBarry Smith static PetscErrorCode PCASMSetType_ASM(PC pc,PCASMType type) 8014b9ad928SBarry Smith { 80292bb6962SLisandro Dalcin PC_ASM *osm = (PC_ASM*)pc->data; 8034b9ad928SBarry Smith 8044b9ad928SBarry Smith PetscFunctionBegin; 8054b9ad928SBarry Smith osm->type = type; 806bf108f30SBarry Smith osm->type_set = PETSC_TRUE; 8074b9ad928SBarry Smith PetscFunctionReturn(0); 8084b9ad928SBarry Smith } 8094b9ad928SBarry Smith 810c60c7ad4SBarry Smith static PetscErrorCode PCASMGetType_ASM(PC pc,PCASMType *type) 811c60c7ad4SBarry Smith { 812c60c7ad4SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 813c60c7ad4SBarry Smith 814c60c7ad4SBarry Smith PetscFunctionBegin; 815c60c7ad4SBarry Smith *type = osm->type; 816c60c7ad4SBarry Smith PetscFunctionReturn(0); 817c60c7ad4SBarry Smith } 818c60c7ad4SBarry Smith 81912cd4985SMatthew G. Knepley static PetscErrorCode PCASMSetLocalType_ASM(PC pc, PCCompositeType type) 82012cd4985SMatthew G. Knepley { 82112cd4985SMatthew G. Knepley PC_ASM *osm = (PC_ASM *) pc->data; 82212cd4985SMatthew G. Knepley 82312cd4985SMatthew G. Knepley PetscFunctionBegin; 824d0ecd4dfSBarry 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"); 82512cd4985SMatthew G. Knepley osm->loctype = type; 82612cd4985SMatthew G. Knepley PetscFunctionReturn(0); 82712cd4985SMatthew G. Knepley } 82812cd4985SMatthew G. Knepley 82912cd4985SMatthew G. Knepley static PetscErrorCode PCASMGetLocalType_ASM(PC pc, PCCompositeType *type) 83012cd4985SMatthew G. Knepley { 83112cd4985SMatthew G. Knepley PC_ASM *osm = (PC_ASM *) pc->data; 83212cd4985SMatthew G. Knepley 83312cd4985SMatthew G. Knepley PetscFunctionBegin; 83412cd4985SMatthew G. Knepley *type = osm->loctype; 83512cd4985SMatthew G. Knepley PetscFunctionReturn(0); 83612cd4985SMatthew G. Knepley } 83712cd4985SMatthew G. Knepley 8381e6b0712SBarry Smith static PetscErrorCode PCASMSetSortIndices_ASM(PC pc,PetscBool doSort) 8396ed231c7SMatthew Knepley { 8406ed231c7SMatthew Knepley PC_ASM *osm = (PC_ASM*)pc->data; 8416ed231c7SMatthew Knepley 8426ed231c7SMatthew Knepley PetscFunctionBegin; 8436ed231c7SMatthew Knepley osm->sort_indices = doSort; 8446ed231c7SMatthew Knepley PetscFunctionReturn(0); 8456ed231c7SMatthew Knepley } 8466ed231c7SMatthew Knepley 8471e6b0712SBarry Smith static PetscErrorCode PCASMGetSubKSP_ASM(PC pc,PetscInt *n_local,PetscInt *first_local,KSP **ksp) 8484b9ad928SBarry Smith { 84992bb6962SLisandro Dalcin PC_ASM *osm = (PC_ASM*)pc->data; 850dfbe8321SBarry Smith PetscErrorCode ierr; 8514b9ad928SBarry Smith 8524b9ad928SBarry Smith PetscFunctionBegin; 853ce94432eSBarry 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"); 8544b9ad928SBarry Smith 8552fa5cd67SKarl Rupp if (n_local) *n_local = osm->n_local_true; 85692bb6962SLisandro Dalcin if (first_local) { 857ce94432eSBarry Smith ierr = MPI_Scan(&osm->n_local_true,first_local,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr); 85892bb6962SLisandro Dalcin *first_local -= osm->n_local_true; 85992bb6962SLisandro Dalcin } 86092bb6962SLisandro Dalcin if (ksp) { 86192bb6962SLisandro Dalcin /* Assume that local solves are now different; not necessarily 86292bb6962SLisandro Dalcin true though! This flag is used only for PCView_ASM() */ 86392bb6962SLisandro Dalcin *ksp = osm->ksp; 86492bb6962SLisandro Dalcin osm->same_local_solves = PETSC_FALSE; 86592bb6962SLisandro Dalcin } 8664b9ad928SBarry Smith PetscFunctionReturn(0); 8674b9ad928SBarry Smith } 8684b9ad928SBarry Smith 86980ec0b7dSPatrick Sanan static PetscErrorCode PCASMGetSubMatType_ASM(PC pc,MatType *sub_mat_type) 87080ec0b7dSPatrick Sanan { 87180ec0b7dSPatrick Sanan PC_ASM *osm = (PC_ASM*)pc->data; 87280ec0b7dSPatrick Sanan 87380ec0b7dSPatrick Sanan PetscFunctionBegin; 87480ec0b7dSPatrick Sanan PetscValidHeaderSpecific(pc,PC_CLASSID,1); 87580ec0b7dSPatrick Sanan PetscValidPointer(sub_mat_type,2); 87680ec0b7dSPatrick Sanan *sub_mat_type = osm->sub_mat_type; 87780ec0b7dSPatrick Sanan PetscFunctionReturn(0); 87880ec0b7dSPatrick Sanan } 87980ec0b7dSPatrick Sanan 88080ec0b7dSPatrick Sanan static PetscErrorCode PCASMSetSubMatType_ASM(PC pc,MatType sub_mat_type) 88180ec0b7dSPatrick Sanan { 88280ec0b7dSPatrick Sanan PetscErrorCode ierr; 88380ec0b7dSPatrick Sanan PC_ASM *osm = (PC_ASM*)pc->data; 88480ec0b7dSPatrick Sanan 88580ec0b7dSPatrick Sanan PetscFunctionBegin; 88680ec0b7dSPatrick Sanan PetscValidHeaderSpecific(pc,PC_CLASSID,1); 88780ec0b7dSPatrick Sanan ierr = PetscFree(osm->sub_mat_type);CHKERRQ(ierr); 88880ec0b7dSPatrick Sanan ierr = PetscStrallocpy(sub_mat_type,(char**)&osm->sub_mat_type);CHKERRQ(ierr); 88980ec0b7dSPatrick Sanan PetscFunctionReturn(0); 89080ec0b7dSPatrick Sanan } 89180ec0b7dSPatrick Sanan 8924b9ad928SBarry Smith /*@C 8931093a601SBarry Smith PCASMSetLocalSubdomains - Sets the local subdomains (for this processor only) for the additive Schwarz preconditioner. 8944b9ad928SBarry Smith 8954b9ad928SBarry Smith Collective on PC 8964b9ad928SBarry Smith 8974b9ad928SBarry Smith Input Parameters: 8984b9ad928SBarry Smith + pc - the preconditioner context 8994b9ad928SBarry Smith . n - the number of subdomains for this processor (default value = 1) 9008c03b21aSDmitry Karpeev . is - the index set that defines the subdomains for this processor 9010298fd71SBarry Smith (or NULL for PETSc to determine subdomains) 902f1ee410cSBarry 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 903f1ee410cSBarry Smith (or NULL to not provide these) 9044b9ad928SBarry Smith 905342c94f9SMatthew G. Knepley Options Database Key: 906342c94f9SMatthew G. Knepley To set the total number of subdomain blocks rather than specify the 907342c94f9SMatthew G. Knepley index sets, use the option 908342c94f9SMatthew G. Knepley . -pc_asm_local_blocks <blks> - Sets local blocks 909342c94f9SMatthew G. Knepley 9104b9ad928SBarry Smith Notes: 9111093a601SBarry Smith The IS numbering is in the parallel, global numbering of the vector for both is and is_local 9124b9ad928SBarry Smith 9134b9ad928SBarry Smith By default the ASM preconditioner uses 1 block per processor. 9144b9ad928SBarry Smith 9154b9ad928SBarry Smith Use PCASMSetTotalSubdomains() to set the subdomains for all processors. 9164b9ad928SBarry Smith 917f1ee410cSBarry Smith If is_local is provided and PCASMType is PC_ASM_RESTRICT then the solution only over the is_local region is interpolated 918f1ee410cSBarry Smith back to form the global solution (this is the standard restricted additive Schwarz method) 919f1ee410cSBarry Smith 920f1ee410cSBarry 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 921f1ee410cSBarry Smith no code to handle that case. 922f1ee410cSBarry Smith 9234b9ad928SBarry Smith Level: advanced 9244b9ad928SBarry Smith 9254b9ad928SBarry Smith .keywords: PC, ASM, set, local, subdomains, additive Schwarz 9264b9ad928SBarry Smith 9274b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(), 928f1ee410cSBarry Smith PCASMCreateSubdomains2D(), PCASMGetLocalSubdomains(), PCASMType, PCASMSetType() 9294b9ad928SBarry Smith @*/ 9307087cfbeSBarry Smith PetscErrorCode PCASMSetLocalSubdomains(PC pc,PetscInt n,IS is[],IS is_local[]) 9314b9ad928SBarry Smith { 9327bb14e67SBarry Smith PetscErrorCode ierr; 9334b9ad928SBarry Smith 9344b9ad928SBarry Smith PetscFunctionBegin; 9350700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 9367bb14e67SBarry Smith ierr = PetscTryMethod(pc,"PCASMSetLocalSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,n,is,is_local));CHKERRQ(ierr); 9374b9ad928SBarry Smith PetscFunctionReturn(0); 9384b9ad928SBarry Smith } 9394b9ad928SBarry Smith 9404b9ad928SBarry Smith /*@C 941feb221f8SDmitry Karpeev PCASMSetTotalSubdomains - Sets the subdomains for all processors for the 9424b9ad928SBarry Smith additive Schwarz preconditioner. Either all or no processors in the 9434b9ad928SBarry Smith PC communicator must call this routine, with the same index sets. 9444b9ad928SBarry Smith 9454b9ad928SBarry Smith Collective on PC 9464b9ad928SBarry Smith 9474b9ad928SBarry Smith Input Parameters: 9484b9ad928SBarry Smith + pc - the preconditioner context 949feb221f8SDmitry Karpeev . N - the number of subdomains for all processors 950feb221f8SDmitry Karpeev . is - the index sets that define the subdomains for all processors 951dfaa0c5fSBarry Smith (or NULL to ask PETSc to determine the subdomains) 9522b691e39SMatthew Knepley - is_local - the index sets that define the local part of the subdomains for this processor 953f1ee410cSBarry Smith (or NULL to not provide this information) 9544b9ad928SBarry Smith 9554b9ad928SBarry Smith Options Database Key: 9564b9ad928SBarry Smith To set the total number of subdomain blocks rather than specify the 9574b9ad928SBarry Smith index sets, use the option 9584b9ad928SBarry Smith . -pc_asm_blocks <blks> - Sets total blocks 9594b9ad928SBarry Smith 9604b9ad928SBarry Smith Notes: 961f1ee410cSBarry Smith Currently you cannot use this to set the actual subdomains with the argument is or is_local. 9624b9ad928SBarry Smith 9634b9ad928SBarry Smith By default the ASM preconditioner uses 1 block per processor. 9644b9ad928SBarry Smith 9654b9ad928SBarry Smith These index sets cannot be destroyed until after completion of the 9664b9ad928SBarry Smith linear solves for which the ASM preconditioner is being used. 9674b9ad928SBarry Smith 9684b9ad928SBarry Smith Use PCASMSetLocalSubdomains() to set local subdomains. 9694b9ad928SBarry Smith 9701093a601SBarry Smith The IS numbering is in the parallel, global numbering of the vector for both is and is_local 9711093a601SBarry Smith 9724b9ad928SBarry Smith Level: advanced 9734b9ad928SBarry Smith 9744b9ad928SBarry Smith .keywords: PC, ASM, set, total, global, subdomains, additive Schwarz 9754b9ad928SBarry Smith 9764b9ad928SBarry Smith .seealso: PCASMSetLocalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(), 9774b9ad928SBarry Smith PCASMCreateSubdomains2D() 9784b9ad928SBarry Smith @*/ 9797087cfbeSBarry Smith PetscErrorCode PCASMSetTotalSubdomains(PC pc,PetscInt N,IS is[],IS is_local[]) 9804b9ad928SBarry Smith { 9817bb14e67SBarry Smith PetscErrorCode ierr; 9824b9ad928SBarry Smith 9834b9ad928SBarry Smith PetscFunctionBegin; 9840700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 9857bb14e67SBarry Smith ierr = PetscTryMethod(pc,"PCASMSetTotalSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,N,is,is_local));CHKERRQ(ierr); 9864b9ad928SBarry Smith PetscFunctionReturn(0); 9874b9ad928SBarry Smith } 9884b9ad928SBarry Smith 9894b9ad928SBarry Smith /*@ 9904b9ad928SBarry Smith PCASMSetOverlap - Sets the overlap between a pair of subdomains for the 9914b9ad928SBarry Smith additive Schwarz preconditioner. Either all or no processors in the 992f1ee410cSBarry Smith PC communicator must call this routine. 9934b9ad928SBarry Smith 994ad4df100SBarry Smith Logically Collective on PC 9954b9ad928SBarry Smith 9964b9ad928SBarry Smith Input Parameters: 9974b9ad928SBarry Smith + pc - the preconditioner context 9984b9ad928SBarry Smith - ovl - the amount of overlap between subdomains (ovl >= 0, default value = 1) 9994b9ad928SBarry Smith 10004b9ad928SBarry Smith Options Database Key: 10014b9ad928SBarry Smith . -pc_asm_overlap <ovl> - Sets overlap 10024b9ad928SBarry Smith 10034b9ad928SBarry Smith Notes: 10044b9ad928SBarry Smith By default the ASM preconditioner uses 1 block per processor. To use 10054b9ad928SBarry Smith multiple blocks per perocessor, see PCASMSetTotalSubdomains() and 10064b9ad928SBarry Smith PCASMSetLocalSubdomains() (and the option -pc_asm_blocks <blks>). 10074b9ad928SBarry Smith 10084b9ad928SBarry Smith The overlap defaults to 1, so if one desires that no additional 10094b9ad928SBarry Smith overlap be computed beyond what may have been set with a call to 10104b9ad928SBarry Smith PCASMSetTotalSubdomains() or PCASMSetLocalSubdomains(), then ovl 10114b9ad928SBarry Smith must be set to be 0. In particular, if one does not explicitly set 10124b9ad928SBarry Smith the subdomains an application code, then all overlap would be computed 10134b9ad928SBarry Smith internally by PETSc, and using an overlap of 0 would result in an ASM 10144b9ad928SBarry Smith variant that is equivalent to the block Jacobi preconditioner. 10154b9ad928SBarry Smith 1016f1ee410cSBarry Smith The default algorithm used by PETSc to increase overlap is fast, but not scalable, 1017f1ee410cSBarry Smith use the option -mat_increase_overlap_scalable when the problem and number of processes is large. 1018f1ee410cSBarry Smith 10194b9ad928SBarry Smith Note that one can define initial index sets with any overlap via 1020f1ee410cSBarry Smith PCASMSetLocalSubdomains(); the routine 10214b9ad928SBarry Smith PCASMSetOverlap() merely allows PETSc to extend that overlap further 10224b9ad928SBarry Smith if desired. 10234b9ad928SBarry Smith 10244b9ad928SBarry Smith Level: intermediate 10254b9ad928SBarry Smith 10264b9ad928SBarry Smith .keywords: PC, ASM, set, overlap 10274b9ad928SBarry Smith 10284b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetLocalSubdomains(), PCASMGetSubKSP(), 1029f1ee410cSBarry Smith PCASMCreateSubdomains2D(), PCASMGetLocalSubdomains(), MatIncreaseOverlap() 10304b9ad928SBarry Smith @*/ 10317087cfbeSBarry Smith PetscErrorCode PCASMSetOverlap(PC pc,PetscInt ovl) 10324b9ad928SBarry Smith { 10337bb14e67SBarry Smith PetscErrorCode ierr; 10344b9ad928SBarry Smith 10354b9ad928SBarry Smith PetscFunctionBegin; 10360700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1037c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,ovl,2); 10387bb14e67SBarry Smith ierr = PetscTryMethod(pc,"PCASMSetOverlap_C",(PC,PetscInt),(pc,ovl));CHKERRQ(ierr); 10394b9ad928SBarry Smith PetscFunctionReturn(0); 10404b9ad928SBarry Smith } 10414b9ad928SBarry Smith 10424b9ad928SBarry Smith /*@ 10434b9ad928SBarry Smith PCASMSetType - Sets the type of restriction and interpolation used 10444b9ad928SBarry Smith for local problems in the additive Schwarz method. 10454b9ad928SBarry Smith 1046ad4df100SBarry Smith Logically Collective on PC 10474b9ad928SBarry Smith 10484b9ad928SBarry Smith Input Parameters: 10494b9ad928SBarry Smith + pc - the preconditioner context 10504b9ad928SBarry Smith - type - variant of ASM, one of 10514b9ad928SBarry Smith .vb 10524b9ad928SBarry Smith PC_ASM_BASIC - full interpolation and restriction 10534b9ad928SBarry Smith PC_ASM_RESTRICT - full restriction, local processor interpolation 10544b9ad928SBarry Smith PC_ASM_INTERPOLATE - full interpolation, local processor restriction 10554b9ad928SBarry Smith PC_ASM_NONE - local processor restriction and interpolation 10564b9ad928SBarry Smith .ve 10574b9ad928SBarry Smith 10584b9ad928SBarry Smith Options Database Key: 10594b9ad928SBarry Smith . -pc_asm_type [basic,restrict,interpolate,none] - Sets ASM type 10604b9ad928SBarry Smith 106195452b02SPatrick Sanan Notes: 106295452b02SPatrick Sanan if the is_local arguments are passed to PCASMSetLocalSubdomains() then they are used when PC_ASM_RESTRICT has been selected 1063f1ee410cSBarry Smith to limit the local processor interpolation 1064f1ee410cSBarry Smith 10654b9ad928SBarry Smith Level: intermediate 10664b9ad928SBarry Smith 10674b9ad928SBarry Smith .keywords: PC, ASM, set, type 10684b9ad928SBarry Smith 10694b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetTotalSubdomains(), PCASMGetSubKSP(), 1070f1ee410cSBarry Smith PCASMCreateSubdomains2D(), PCASMType, PCASMSetLocalType(), PCASMGetLocalType() 10714b9ad928SBarry Smith @*/ 10727087cfbeSBarry Smith PetscErrorCode PCASMSetType(PC pc,PCASMType type) 10734b9ad928SBarry Smith { 10747bb14e67SBarry Smith PetscErrorCode ierr; 10754b9ad928SBarry Smith 10764b9ad928SBarry Smith PetscFunctionBegin; 10770700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1078c5eb9154SBarry Smith PetscValidLogicalCollectiveEnum(pc,type,2); 10797bb14e67SBarry Smith ierr = PetscTryMethod(pc,"PCASMSetType_C",(PC,PCASMType),(pc,type));CHKERRQ(ierr); 10804b9ad928SBarry Smith PetscFunctionReturn(0); 10814b9ad928SBarry Smith } 10824b9ad928SBarry Smith 1083c60c7ad4SBarry Smith /*@ 1084c60c7ad4SBarry Smith PCASMGetType - Gets the type of restriction and interpolation used 1085c60c7ad4SBarry Smith for local problems in the additive Schwarz method. 1086c60c7ad4SBarry Smith 1087c60c7ad4SBarry Smith Logically Collective on PC 1088c60c7ad4SBarry Smith 1089c60c7ad4SBarry Smith Input Parameter: 1090c60c7ad4SBarry Smith . pc - the preconditioner context 1091c60c7ad4SBarry Smith 1092c60c7ad4SBarry Smith Output Parameter: 1093c60c7ad4SBarry Smith . type - variant of ASM, one of 1094c60c7ad4SBarry Smith 1095c60c7ad4SBarry Smith .vb 1096c60c7ad4SBarry Smith PC_ASM_BASIC - full interpolation and restriction 1097c60c7ad4SBarry Smith PC_ASM_RESTRICT - full restriction, local processor interpolation 1098c60c7ad4SBarry Smith PC_ASM_INTERPOLATE - full interpolation, local processor restriction 1099c60c7ad4SBarry Smith PC_ASM_NONE - local processor restriction and interpolation 1100c60c7ad4SBarry Smith .ve 1101c60c7ad4SBarry Smith 1102c60c7ad4SBarry Smith Options Database Key: 1103c60c7ad4SBarry Smith . -pc_asm_type [basic,restrict,interpolate,none] - Sets ASM type 1104c60c7ad4SBarry Smith 1105c60c7ad4SBarry Smith Level: intermediate 1106c60c7ad4SBarry Smith 1107c60c7ad4SBarry Smith .keywords: PC, ASM, set, type 1108c60c7ad4SBarry Smith 1109c60c7ad4SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetTotalSubdomains(), PCASMGetSubKSP(), 1110f1ee410cSBarry Smith PCASMCreateSubdomains2D(), PCASMType, PCASMSetType(), PCASMSetLocalType(), PCASMGetLocalType() 1111c60c7ad4SBarry Smith @*/ 1112c60c7ad4SBarry Smith PetscErrorCode PCASMGetType(PC pc,PCASMType *type) 1113c60c7ad4SBarry Smith { 1114c60c7ad4SBarry Smith PetscErrorCode ierr; 1115c60c7ad4SBarry Smith 1116c60c7ad4SBarry Smith PetscFunctionBegin; 1117c60c7ad4SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1118c60c7ad4SBarry Smith ierr = PetscUseMethod(pc,"PCASMGetType_C",(PC,PCASMType*),(pc,type));CHKERRQ(ierr); 1119c60c7ad4SBarry Smith PetscFunctionReturn(0); 1120c60c7ad4SBarry Smith } 1121c60c7ad4SBarry Smith 112212cd4985SMatthew G. Knepley /*@ 112312cd4985SMatthew G. Knepley PCASMSetLocalType - Sets the type of composition used for local problems in the additive Schwarz method. 112412cd4985SMatthew G. Knepley 112512cd4985SMatthew G. Knepley Logically Collective on PC 112612cd4985SMatthew G. Knepley 112712cd4985SMatthew G. Knepley Input Parameters: 112812cd4985SMatthew G. Knepley + pc - the preconditioner context 112912cd4985SMatthew G. Knepley - type - type of composition, one of 113012cd4985SMatthew G. Knepley .vb 113112cd4985SMatthew G. Knepley PC_COMPOSITE_ADDITIVE - local additive combination 113212cd4985SMatthew G. Knepley PC_COMPOSITE_MULTIPLICATIVE - local multiplicative combination 113312cd4985SMatthew G. Knepley .ve 113412cd4985SMatthew G. Knepley 113512cd4985SMatthew G. Knepley Options Database Key: 113612cd4985SMatthew G. Knepley . -pc_asm_local_type [additive,multiplicative] - Sets local solver composition type 113712cd4985SMatthew G. Knepley 113812cd4985SMatthew G. Knepley Level: intermediate 113912cd4985SMatthew G. Knepley 1140f1ee410cSBarry Smith .seealso: PCASMSetType(), PCASMGetType(), PCASMGetLocalType(), PCASM, PCASMType, PCASMSetType(), PCASMGetType(), PCCompositeType 114112cd4985SMatthew G. Knepley @*/ 114212cd4985SMatthew G. Knepley PetscErrorCode PCASMSetLocalType(PC pc, PCCompositeType type) 114312cd4985SMatthew G. Knepley { 114412cd4985SMatthew G. Knepley PetscErrorCode ierr; 114512cd4985SMatthew G. Knepley 114612cd4985SMatthew G. Knepley PetscFunctionBegin; 114712cd4985SMatthew G. Knepley PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 114812cd4985SMatthew G. Knepley PetscValidLogicalCollectiveEnum(pc, type, 2); 114912cd4985SMatthew G. Knepley ierr = PetscTryMethod(pc, "PCASMSetLocalType_C", (PC, PCCompositeType), (pc, type));CHKERRQ(ierr); 115012cd4985SMatthew G. Knepley PetscFunctionReturn(0); 115112cd4985SMatthew G. Knepley } 115212cd4985SMatthew G. Knepley 115312cd4985SMatthew G. Knepley /*@ 115412cd4985SMatthew G. Knepley PCASMGetLocalType - Gets the type of composition used for local problems in the additive Schwarz method. 115512cd4985SMatthew G. Knepley 115612cd4985SMatthew G. Knepley Logically Collective on PC 115712cd4985SMatthew G. Knepley 115812cd4985SMatthew G. Knepley Input Parameter: 115912cd4985SMatthew G. Knepley . pc - the preconditioner context 116012cd4985SMatthew G. Knepley 116112cd4985SMatthew G. Knepley Output Parameter: 116212cd4985SMatthew G. Knepley . type - type of composition, one of 116312cd4985SMatthew G. Knepley .vb 116412cd4985SMatthew G. Knepley PC_COMPOSITE_ADDITIVE - local additive combination 116512cd4985SMatthew G. Knepley PC_COMPOSITE_MULTIPLICATIVE - local multiplicative combination 116612cd4985SMatthew G. Knepley .ve 116712cd4985SMatthew G. Knepley 116812cd4985SMatthew G. Knepley Options Database Key: 116912cd4985SMatthew G. Knepley . -pc_asm_local_type [additive,multiplicative] - Sets local solver composition type 117012cd4985SMatthew G. Knepley 117112cd4985SMatthew G. Knepley Level: intermediate 117212cd4985SMatthew G. Knepley 1173f1ee410cSBarry Smith .seealso: PCASMSetType(), PCASMGetType(), PCASMSetLocalType(), PCASMCreate(), PCASMType, PCASMSetType(), PCASMGetType(), PCCompositeType 117412cd4985SMatthew G. Knepley @*/ 117512cd4985SMatthew G. Knepley PetscErrorCode PCASMGetLocalType(PC pc, PCCompositeType *type) 117612cd4985SMatthew G. Knepley { 117712cd4985SMatthew G. Knepley PetscErrorCode ierr; 117812cd4985SMatthew G. Knepley 117912cd4985SMatthew G. Knepley PetscFunctionBegin; 118012cd4985SMatthew G. Knepley PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 118112cd4985SMatthew G. Knepley PetscValidPointer(type, 2); 118212cd4985SMatthew G. Knepley ierr = PetscUseMethod(pc, "PCASMGetLocalType_C", (PC, PCCompositeType *), (pc, type));CHKERRQ(ierr); 118312cd4985SMatthew G. Knepley PetscFunctionReturn(0); 118412cd4985SMatthew G. Knepley } 118512cd4985SMatthew G. Knepley 11866ed231c7SMatthew Knepley /*@ 11876ed231c7SMatthew Knepley PCASMSetSortIndices - Determines whether subdomain indices are sorted. 11886ed231c7SMatthew Knepley 1189ad4df100SBarry Smith Logically Collective on PC 11906ed231c7SMatthew Knepley 11916ed231c7SMatthew Knepley Input Parameters: 11926ed231c7SMatthew Knepley + pc - the preconditioner context 11936ed231c7SMatthew Knepley - doSort - sort the subdomain indices 11946ed231c7SMatthew Knepley 11956ed231c7SMatthew Knepley Level: intermediate 11966ed231c7SMatthew Knepley 11976ed231c7SMatthew Knepley .keywords: PC, ASM, set, type 11986ed231c7SMatthew Knepley 11996ed231c7SMatthew Knepley .seealso: PCASMSetLocalSubdomains(), PCASMSetTotalSubdomains(), PCASMGetSubKSP(), 12006ed231c7SMatthew Knepley PCASMCreateSubdomains2D() 12016ed231c7SMatthew Knepley @*/ 12027087cfbeSBarry Smith PetscErrorCode PCASMSetSortIndices(PC pc,PetscBool doSort) 12036ed231c7SMatthew Knepley { 12047bb14e67SBarry Smith PetscErrorCode ierr; 12056ed231c7SMatthew Knepley 12066ed231c7SMatthew Knepley PetscFunctionBegin; 12070700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1208acfcf0e5SJed Brown PetscValidLogicalCollectiveBool(pc,doSort,2); 12097bb14e67SBarry Smith ierr = PetscTryMethod(pc,"PCASMSetSortIndices_C",(PC,PetscBool),(pc,doSort));CHKERRQ(ierr); 12106ed231c7SMatthew Knepley PetscFunctionReturn(0); 12116ed231c7SMatthew Knepley } 12126ed231c7SMatthew Knepley 12134b9ad928SBarry Smith /*@C 12144b9ad928SBarry Smith PCASMGetSubKSP - Gets the local KSP contexts for all blocks on 12154b9ad928SBarry Smith this processor. 12164b9ad928SBarry Smith 12174b9ad928SBarry Smith Collective on PC iff first_local is requested 12184b9ad928SBarry Smith 12194b9ad928SBarry Smith Input Parameter: 12204b9ad928SBarry Smith . pc - the preconditioner context 12214b9ad928SBarry Smith 12224b9ad928SBarry Smith Output Parameters: 12230298fd71SBarry Smith + n_local - the number of blocks on this processor or NULL 12240298fd71SBarry Smith . first_local - the global number of the first block on this processor or NULL, 12250298fd71SBarry Smith all processors must request or all must pass NULL 12264b9ad928SBarry Smith - ksp - the array of KSP contexts 12274b9ad928SBarry Smith 12284b9ad928SBarry Smith Note: 1229d29017ddSJed Brown After PCASMGetSubKSP() the array of KSPes is not to be freed. 12304b9ad928SBarry Smith 12314b9ad928SBarry Smith You must call KSPSetUp() before calling PCASMGetSubKSP(). 12324b9ad928SBarry Smith 1233d29017ddSJed Brown Fortran note: 12342bf68e3eSBarry 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. 1235d29017ddSJed Brown 12364b9ad928SBarry Smith Level: advanced 12374b9ad928SBarry Smith 12384b9ad928SBarry Smith .keywords: PC, ASM, additive Schwarz, get, sub, KSP, context 12394b9ad928SBarry Smith 12404b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetTotalSubdomains(), PCASMSetOverlap(), 12414b9ad928SBarry Smith PCASMCreateSubdomains2D(), 12424b9ad928SBarry Smith @*/ 12437087cfbeSBarry Smith PetscErrorCode PCASMGetSubKSP(PC pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp[]) 12444b9ad928SBarry Smith { 12457bb14e67SBarry Smith PetscErrorCode ierr; 12464b9ad928SBarry Smith 12474b9ad928SBarry Smith PetscFunctionBegin; 12480700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 12497bb14e67SBarry Smith ierr = PetscUseMethod(pc,"PCASMGetSubKSP_C",(PC,PetscInt*,PetscInt*,KSP **),(pc,n_local,first_local,ksp));CHKERRQ(ierr); 12504b9ad928SBarry Smith PetscFunctionReturn(0); 12514b9ad928SBarry Smith } 12524b9ad928SBarry Smith 12534b9ad928SBarry Smith /* -------------------------------------------------------------------------------------*/ 12544b9ad928SBarry Smith /*MC 12553b09bd56SBarry Smith PCASM - Use the (restricted) additive Schwarz method, each block is (approximately) solved with 12564b9ad928SBarry Smith its own KSP object. 12574b9ad928SBarry Smith 12584b9ad928SBarry Smith Options Database Keys: 125949517cdeSBarry Smith + -pc_asm_blocks <blks> - Sets total blocks 12604b9ad928SBarry Smith . -pc_asm_overlap <ovl> - Sets overlap 1261f1ee410cSBarry Smith . -pc_asm_type [basic,restrict,interpolate,none] - Sets ASM type, default is restrict 1262f1ee410cSBarry Smith - -pc_asm_local_type [additive, multiplicative] - Sets ASM type, default is additive 12634b9ad928SBarry Smith 12643b09bd56SBarry Smith IMPORTANT: If you run with, for example, 3 blocks on 1 processor or 3 blocks on 3 processors you 12653b09bd56SBarry Smith will get a different convergence rate due to the default option of -pc_asm_type restrict. Use 12663b09bd56SBarry Smith -pc_asm_type basic to use the standard ASM. 12673b09bd56SBarry Smith 126895452b02SPatrick Sanan Notes: 126995452b02SPatrick Sanan Each processor can have one or more blocks, but a block cannot be shared by more 1270f1ee410cSBarry Smith than one processor. Use PCGASM for subdomains shared by multiple processes. Defaults to one block per processor. 12714b9ad928SBarry Smith 12723b09bd56SBarry Smith To set options on the solvers for each block append -sub_ to all the KSP, and PC 1273d7ee0231SBarry Smith options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly 12744b9ad928SBarry Smith 1275a8c7a070SBarry Smith To set the options on the solvers separate for each block call PCASMGetSubKSP() 12764b9ad928SBarry Smith and set the options directly on the resulting KSP object (you can access its PC 12774b9ad928SBarry Smith with KSPGetPC()) 12784b9ad928SBarry Smith 12794b9ad928SBarry Smith Level: beginner 12804b9ad928SBarry Smith 12814b9ad928SBarry Smith Concepts: additive Schwarz method 12824b9ad928SBarry Smith 1283c582cd25SBarry Smith References: 128496a0c994SBarry Smith + 1. - M Dryja, OB Widlund, An additive variant of the Schwarz alternating method for the case of many subregions 128596a0c994SBarry Smith Courant Institute, New York University Technical report 128696a0c994SBarry Smith - 1. - Barry Smith, Petter Bjorstad, and William Gropp, Domain Decompositions: Parallel Multilevel Methods for Elliptic Partial Differential Equations, 128796a0c994SBarry Smith Cambridge University Press. 1288c582cd25SBarry Smith 12894b9ad928SBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 1290f1ee410cSBarry Smith PCBJACOBI, PCASMGetSubKSP(), PCASMSetLocalSubdomains(), PCASMType, PCASMGetType(), PCASMSetLocalType(), PCASMGetLocalType() 1291f1ee410cSBarry Smith PCASMSetTotalSubdomains(), PCSetModifySubmatrices(), PCASMSetOverlap(), PCASMSetType(), PCCompositeType 1292e09e08ccSBarry Smith 12934b9ad928SBarry Smith M*/ 12944b9ad928SBarry Smith 12958cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_ASM(PC pc) 12964b9ad928SBarry Smith { 1297dfbe8321SBarry Smith PetscErrorCode ierr; 12984b9ad928SBarry Smith PC_ASM *osm; 12994b9ad928SBarry Smith 13004b9ad928SBarry Smith PetscFunctionBegin; 1301b00a9115SJed Brown ierr = PetscNewLog(pc,&osm);CHKERRQ(ierr); 13022fa5cd67SKarl Rupp 13034b9ad928SBarry Smith osm->n = PETSC_DECIDE; 13044b9ad928SBarry Smith osm->n_local = 0; 13052b837212SDmitry Karpeev osm->n_local_true = PETSC_DECIDE; 13064b9ad928SBarry Smith osm->overlap = 1; 13074b9ad928SBarry Smith osm->ksp = 0; 13082b691e39SMatthew Knepley osm->restriction = 0; 13095b3c0d42Seaulisa osm->lprolongation = 0; 13101dd8081eSeaulisa osm->lrestriction = 0; 131192bb6962SLisandro Dalcin osm->x = 0; 131292bb6962SLisandro Dalcin osm->y = 0; 13134b9ad928SBarry Smith osm->is = 0; 13142b691e39SMatthew Knepley osm->is_local = 0; 13154b9ad928SBarry Smith osm->mat = 0; 13164b9ad928SBarry Smith osm->pmat = 0; 13174b9ad928SBarry Smith osm->type = PC_ASM_RESTRICT; 131812cd4985SMatthew G. Knepley osm->loctype = PC_COMPOSITE_ADDITIVE; 13194b9ad928SBarry Smith osm->same_local_solves = PETSC_TRUE; 13206ed231c7SMatthew Knepley osm->sort_indices = PETSC_TRUE; 1321d709ea83SDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 132280ec0b7dSPatrick Sanan osm->sub_mat_type = NULL; 13234b9ad928SBarry Smith 132492bb6962SLisandro Dalcin pc->data = (void*)osm; 13254b9ad928SBarry Smith pc->ops->apply = PCApply_ASM; 13264b9ad928SBarry Smith pc->ops->applytranspose = PCApplyTranspose_ASM; 13274b9ad928SBarry Smith pc->ops->setup = PCSetUp_ASM; 1328e91c6855SBarry Smith pc->ops->reset = PCReset_ASM; 13294b9ad928SBarry Smith pc->ops->destroy = PCDestroy_ASM; 13304b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_ASM; 13314b9ad928SBarry Smith pc->ops->setuponblocks = PCSetUpOnBlocks_ASM; 13324b9ad928SBarry Smith pc->ops->view = PCView_ASM; 13334b9ad928SBarry Smith pc->ops->applyrichardson = 0; 13344b9ad928SBarry Smith 1335bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetLocalSubdomains_C",PCASMSetLocalSubdomains_ASM);CHKERRQ(ierr); 1336bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetTotalSubdomains_C",PCASMSetTotalSubdomains_ASM);CHKERRQ(ierr); 1337bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetOverlap_C",PCASMSetOverlap_ASM);CHKERRQ(ierr); 1338bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetType_C",PCASMSetType_ASM);CHKERRQ(ierr); 1339c60c7ad4SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMGetType_C",PCASMGetType_ASM);CHKERRQ(ierr); 134012cd4985SMatthew G. Knepley ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetLocalType_C",PCASMSetLocalType_ASM);CHKERRQ(ierr); 134112cd4985SMatthew G. Knepley ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMGetLocalType_C",PCASMGetLocalType_ASM);CHKERRQ(ierr); 1342bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetSortIndices_C",PCASMSetSortIndices_ASM);CHKERRQ(ierr); 1343bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMGetSubKSP_C",PCASMGetSubKSP_ASM);CHKERRQ(ierr); 134480ec0b7dSPatrick Sanan ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMGetSubMatType_C",PCASMGetSubMatType_ASM);CHKERRQ(ierr); 134580ec0b7dSPatrick Sanan ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetSubMatType_C",PCASMSetSubMatType_ASM);CHKERRQ(ierr); 13464b9ad928SBarry Smith PetscFunctionReturn(0); 13474b9ad928SBarry Smith } 13484b9ad928SBarry Smith 134992bb6962SLisandro Dalcin /*@C 135092bb6962SLisandro Dalcin PCASMCreateSubdomains - Creates the index sets for the overlapping Schwarz 135192bb6962SLisandro Dalcin preconditioner for a any problem on a general grid. 135292bb6962SLisandro Dalcin 135392bb6962SLisandro Dalcin Collective 135492bb6962SLisandro Dalcin 135592bb6962SLisandro Dalcin Input Parameters: 135692bb6962SLisandro Dalcin + A - The global matrix operator 135792bb6962SLisandro Dalcin - n - the number of local blocks 135892bb6962SLisandro Dalcin 135992bb6962SLisandro Dalcin Output Parameters: 136092bb6962SLisandro Dalcin . outis - the array of index sets defining the subdomains 136192bb6962SLisandro Dalcin 136292bb6962SLisandro Dalcin Level: advanced 136392bb6962SLisandro Dalcin 13647d6bfa3bSBarry Smith Note: this generates nonoverlapping subdomains; the PCASM will generate the overlap 13657d6bfa3bSBarry Smith from these if you use PCASMSetLocalSubdomains() 13667d6bfa3bSBarry Smith 13677d6bfa3bSBarry Smith In the Fortran version you must provide the array outis[] already allocated of length n. 13687d6bfa3bSBarry Smith 136992bb6962SLisandro Dalcin .keywords: PC, ASM, additive Schwarz, create, subdomains, unstructured grid 137092bb6962SLisandro Dalcin 137192bb6962SLisandro Dalcin .seealso: PCASMSetLocalSubdomains(), PCASMDestroySubdomains() 137292bb6962SLisandro Dalcin @*/ 13737087cfbeSBarry Smith PetscErrorCode PCASMCreateSubdomains(Mat A, PetscInt n, IS* outis[]) 137492bb6962SLisandro Dalcin { 137592bb6962SLisandro Dalcin MatPartitioning mpart; 137692bb6962SLisandro Dalcin const char *prefix; 137792bb6962SLisandro Dalcin PetscInt i,j,rstart,rend,bs; 1378976e8c5aSLisandro Dalcin PetscBool hasop, isbaij = PETSC_FALSE,foundpart = PETSC_FALSE; 13790298fd71SBarry Smith Mat Ad = NULL, adj; 138092bb6962SLisandro Dalcin IS ispart,isnumb,*is; 138192bb6962SLisandro Dalcin PetscErrorCode ierr; 138292bb6962SLisandro Dalcin 138392bb6962SLisandro Dalcin PetscFunctionBegin; 13840700a824SBarry Smith PetscValidHeaderSpecific(A,MAT_CLASSID,1); 138592bb6962SLisandro Dalcin PetscValidPointer(outis,3); 1386c1235816SBarry Smith if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of local blocks must be > 0, n = %D",n); 138792bb6962SLisandro Dalcin 138892bb6962SLisandro Dalcin /* Get prefix, row distribution, and block size */ 138992bb6962SLisandro Dalcin ierr = MatGetOptionsPrefix(A,&prefix);CHKERRQ(ierr); 139092bb6962SLisandro Dalcin ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 139192bb6962SLisandro Dalcin ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr); 139265e19b50SBarry 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); 139365e19b50SBarry Smith 139492bb6962SLisandro Dalcin /* Get diagonal block from matrix if possible */ 1395976e8c5aSLisandro Dalcin ierr = MatHasOperation(A,MATOP_GET_DIAGONAL_BLOCK,&hasop);CHKERRQ(ierr); 1396976e8c5aSLisandro Dalcin if (hasop) { 139711bd1e4dSLisandro Dalcin ierr = MatGetDiagonalBlock(A,&Ad);CHKERRQ(ierr); 139892bb6962SLisandro Dalcin } 139992bb6962SLisandro Dalcin if (Ad) { 1400251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)Ad,MATSEQBAIJ,&isbaij);CHKERRQ(ierr); 1401251f4c67SDmitry Karpeev if (!isbaij) {ierr = PetscObjectTypeCompare((PetscObject)Ad,MATSEQSBAIJ,&isbaij);CHKERRQ(ierr);} 140292bb6962SLisandro Dalcin } 140392bb6962SLisandro Dalcin if (Ad && n > 1) { 1404ace3abfcSBarry Smith PetscBool match,done; 140592bb6962SLisandro Dalcin /* Try to setup a good matrix partitioning if available */ 140692bb6962SLisandro Dalcin ierr = MatPartitioningCreate(PETSC_COMM_SELF,&mpart);CHKERRQ(ierr); 140792bb6962SLisandro Dalcin ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr); 140892bb6962SLisandro Dalcin ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr); 1409251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGCURRENT,&match);CHKERRQ(ierr); 141092bb6962SLisandro Dalcin if (!match) { 1411251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGSQUARE,&match);CHKERRQ(ierr); 141292bb6962SLisandro Dalcin } 141392bb6962SLisandro Dalcin if (!match) { /* assume a "good" partitioner is available */ 14141a83f524SJed Brown PetscInt na; 14151a83f524SJed Brown const PetscInt *ia,*ja; 141692bb6962SLisandro Dalcin ierr = MatGetRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr); 141792bb6962SLisandro Dalcin if (done) { 141892bb6962SLisandro Dalcin /* Build adjacency matrix by hand. Unfortunately a call to 141992bb6962SLisandro Dalcin MatConvert(Ad,MATMPIADJ,MAT_INITIAL_MATRIX,&adj) will 142092bb6962SLisandro Dalcin remove the block-aij structure and we cannot expect 142192bb6962SLisandro Dalcin MatPartitioning to split vertices as we need */ 14221a83f524SJed Brown PetscInt i,j,len,nnz,cnt,*iia=0,*jja=0; 14231a83f524SJed Brown const PetscInt *row; 142492bb6962SLisandro Dalcin nnz = 0; 142592bb6962SLisandro Dalcin for (i=0; i<na; i++) { /* count number of nonzeros */ 142692bb6962SLisandro Dalcin len = ia[i+1] - ia[i]; 142792bb6962SLisandro Dalcin row = ja + ia[i]; 142892bb6962SLisandro Dalcin for (j=0; j<len; j++) { 142992bb6962SLisandro Dalcin if (row[j] == i) { /* don't count diagonal */ 143092bb6962SLisandro Dalcin len--; break; 143192bb6962SLisandro Dalcin } 143292bb6962SLisandro Dalcin } 143392bb6962SLisandro Dalcin nnz += len; 143492bb6962SLisandro Dalcin } 1435854ce69bSBarry Smith ierr = PetscMalloc1(na+1,&iia);CHKERRQ(ierr); 1436854ce69bSBarry Smith ierr = PetscMalloc1(nnz,&jja);CHKERRQ(ierr); 143792bb6962SLisandro Dalcin nnz = 0; 143892bb6962SLisandro Dalcin iia[0] = 0; 143992bb6962SLisandro Dalcin for (i=0; i<na; i++) { /* fill adjacency */ 144092bb6962SLisandro Dalcin cnt = 0; 144192bb6962SLisandro Dalcin len = ia[i+1] - ia[i]; 144292bb6962SLisandro Dalcin row = ja + ia[i]; 144392bb6962SLisandro Dalcin for (j=0; j<len; j++) { 144492bb6962SLisandro Dalcin if (row[j] != i) { /* if not diagonal */ 144592bb6962SLisandro Dalcin jja[nnz+cnt++] = row[j]; 144692bb6962SLisandro Dalcin } 144792bb6962SLisandro Dalcin } 144892bb6962SLisandro Dalcin nnz += cnt; 144992bb6962SLisandro Dalcin iia[i+1] = nnz; 145092bb6962SLisandro Dalcin } 145192bb6962SLisandro Dalcin /* Partitioning of the adjacency matrix */ 14520298fd71SBarry Smith ierr = MatCreateMPIAdj(PETSC_COMM_SELF,na,na,iia,jja,NULL,&adj);CHKERRQ(ierr); 145392bb6962SLisandro Dalcin ierr = MatPartitioningSetAdjacency(mpart,adj);CHKERRQ(ierr); 145492bb6962SLisandro Dalcin ierr = MatPartitioningSetNParts(mpart,n);CHKERRQ(ierr); 145592bb6962SLisandro Dalcin ierr = MatPartitioningApply(mpart,&ispart);CHKERRQ(ierr); 145692bb6962SLisandro Dalcin ierr = ISPartitioningToNumbering(ispart,&isnumb);CHKERRQ(ierr); 1457fcfd50ebSBarry Smith ierr = MatDestroy(&adj);CHKERRQ(ierr); 145892bb6962SLisandro Dalcin foundpart = PETSC_TRUE; 145992bb6962SLisandro Dalcin } 146092bb6962SLisandro Dalcin ierr = MatRestoreRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr); 146192bb6962SLisandro Dalcin } 1462fcfd50ebSBarry Smith ierr = MatPartitioningDestroy(&mpart);CHKERRQ(ierr); 146392bb6962SLisandro Dalcin } 146492bb6962SLisandro Dalcin 1465785e854fSJed Brown ierr = PetscMalloc1(n,&is);CHKERRQ(ierr); 146692bb6962SLisandro Dalcin *outis = is; 146792bb6962SLisandro Dalcin 146892bb6962SLisandro Dalcin if (!foundpart) { 146992bb6962SLisandro Dalcin 147092bb6962SLisandro Dalcin /* Partitioning by contiguous chunks of rows */ 147192bb6962SLisandro Dalcin 147292bb6962SLisandro Dalcin PetscInt mbs = (rend-rstart)/bs; 147392bb6962SLisandro Dalcin PetscInt start = rstart; 147492bb6962SLisandro Dalcin for (i=0; i<n; i++) { 147592bb6962SLisandro Dalcin PetscInt count = (mbs/n + ((mbs % n) > i)) * bs; 147692bb6962SLisandro Dalcin ierr = ISCreateStride(PETSC_COMM_SELF,count,start,1,&is[i]);CHKERRQ(ierr); 147792bb6962SLisandro Dalcin start += count; 147892bb6962SLisandro Dalcin } 147992bb6962SLisandro Dalcin 148092bb6962SLisandro Dalcin } else { 148192bb6962SLisandro Dalcin 148292bb6962SLisandro Dalcin /* Partitioning by adjacency of diagonal block */ 148392bb6962SLisandro Dalcin 148492bb6962SLisandro Dalcin const PetscInt *numbering; 148592bb6962SLisandro Dalcin PetscInt *count,nidx,*indices,*newidx,start=0; 148692bb6962SLisandro Dalcin /* Get node count in each partition */ 1487785e854fSJed Brown ierr = PetscMalloc1(n,&count);CHKERRQ(ierr); 148892bb6962SLisandro Dalcin ierr = ISPartitioningCount(ispart,n,count);CHKERRQ(ierr); 148992bb6962SLisandro Dalcin if (isbaij && bs > 1) { /* adjust for the block-aij case */ 149092bb6962SLisandro Dalcin for (i=0; i<n; i++) count[i] *= bs; 149192bb6962SLisandro Dalcin } 149292bb6962SLisandro Dalcin /* Build indices from node numbering */ 149392bb6962SLisandro Dalcin ierr = ISGetLocalSize(isnumb,&nidx);CHKERRQ(ierr); 1494785e854fSJed Brown ierr = PetscMalloc1(nidx,&indices);CHKERRQ(ierr); 149592bb6962SLisandro Dalcin for (i=0; i<nidx; i++) indices[i] = i; /* needs to be initialized */ 149692bb6962SLisandro Dalcin ierr = ISGetIndices(isnumb,&numbering);CHKERRQ(ierr); 149792bb6962SLisandro Dalcin ierr = PetscSortIntWithPermutation(nidx,numbering,indices);CHKERRQ(ierr); 149892bb6962SLisandro Dalcin ierr = ISRestoreIndices(isnumb,&numbering);CHKERRQ(ierr); 149992bb6962SLisandro Dalcin if (isbaij && bs > 1) { /* adjust for the block-aij case */ 1500785e854fSJed Brown ierr = PetscMalloc1(nidx*bs,&newidx);CHKERRQ(ierr); 15012fa5cd67SKarl Rupp for (i=0; i<nidx; i++) { 15022fa5cd67SKarl Rupp for (j=0; j<bs; j++) newidx[i*bs+j] = indices[i]*bs + j; 15032fa5cd67SKarl Rupp } 150492bb6962SLisandro Dalcin ierr = PetscFree(indices);CHKERRQ(ierr); 150592bb6962SLisandro Dalcin nidx *= bs; 150692bb6962SLisandro Dalcin indices = newidx; 150792bb6962SLisandro Dalcin } 150892bb6962SLisandro Dalcin /* Shift to get global indices */ 150992bb6962SLisandro Dalcin for (i=0; i<nidx; i++) indices[i] += rstart; 151092bb6962SLisandro Dalcin 151192bb6962SLisandro Dalcin /* Build the index sets for each block */ 151292bb6962SLisandro Dalcin for (i=0; i<n; i++) { 151370b3c8c7SBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,count[i],&indices[start],PETSC_COPY_VALUES,&is[i]);CHKERRQ(ierr); 151492bb6962SLisandro Dalcin ierr = ISSort(is[i]);CHKERRQ(ierr); 151592bb6962SLisandro Dalcin start += count[i]; 151692bb6962SLisandro Dalcin } 151792bb6962SLisandro Dalcin 15183bf036e2SBarry Smith ierr = PetscFree(count);CHKERRQ(ierr); 15193bf036e2SBarry Smith ierr = PetscFree(indices);CHKERRQ(ierr); 1520fcfd50ebSBarry Smith ierr = ISDestroy(&isnumb);CHKERRQ(ierr); 1521fcfd50ebSBarry Smith ierr = ISDestroy(&ispart);CHKERRQ(ierr); 152292bb6962SLisandro Dalcin 152392bb6962SLisandro Dalcin } 152492bb6962SLisandro Dalcin PetscFunctionReturn(0); 152592bb6962SLisandro Dalcin } 152692bb6962SLisandro Dalcin 152792bb6962SLisandro Dalcin /*@C 152892bb6962SLisandro Dalcin PCASMDestroySubdomains - Destroys the index sets created with 152992bb6962SLisandro Dalcin PCASMCreateSubdomains(). Should be called after setting subdomains 153092bb6962SLisandro Dalcin with PCASMSetLocalSubdomains(). 153192bb6962SLisandro Dalcin 153292bb6962SLisandro Dalcin Collective 153392bb6962SLisandro Dalcin 153492bb6962SLisandro Dalcin Input Parameters: 153592bb6962SLisandro Dalcin + n - the number of index sets 15362b691e39SMatthew Knepley . is - the array of index sets 15370298fd71SBarry Smith - is_local - the array of local index sets, can be NULL 153892bb6962SLisandro Dalcin 153992bb6962SLisandro Dalcin Level: advanced 154092bb6962SLisandro Dalcin 154192bb6962SLisandro Dalcin .keywords: PC, ASM, additive Schwarz, create, subdomains, unstructured grid 154292bb6962SLisandro Dalcin 154392bb6962SLisandro Dalcin .seealso: PCASMCreateSubdomains(), PCASMSetLocalSubdomains() 154492bb6962SLisandro Dalcin @*/ 15457087cfbeSBarry Smith PetscErrorCode PCASMDestroySubdomains(PetscInt n, IS is[], IS is_local[]) 154692bb6962SLisandro Dalcin { 154792bb6962SLisandro Dalcin PetscInt i; 154892bb6962SLisandro Dalcin PetscErrorCode ierr; 15495fd66863SKarl Rupp 155092bb6962SLisandro Dalcin PetscFunctionBegin; 1551a35b7b57SDmitry Karpeev if (n <= 0) PetscFunctionReturn(0); 1552a35b7b57SDmitry Karpeev if (is) { 155392bb6962SLisandro Dalcin PetscValidPointer(is,2); 1554fcfd50ebSBarry Smith for (i=0; i<n; i++) { ierr = ISDestroy(&is[i]);CHKERRQ(ierr); } 155592bb6962SLisandro Dalcin ierr = PetscFree(is);CHKERRQ(ierr); 1556a35b7b57SDmitry Karpeev } 15572b691e39SMatthew Knepley if (is_local) { 15582b691e39SMatthew Knepley PetscValidPointer(is_local,3); 1559fcfd50ebSBarry Smith for (i=0; i<n; i++) { ierr = ISDestroy(&is_local[i]);CHKERRQ(ierr); } 15602b691e39SMatthew Knepley ierr = PetscFree(is_local);CHKERRQ(ierr); 15612b691e39SMatthew Knepley } 156292bb6962SLisandro Dalcin PetscFunctionReturn(0); 156392bb6962SLisandro Dalcin } 156492bb6962SLisandro Dalcin 15654b9ad928SBarry Smith /*@ 15664b9ad928SBarry Smith PCASMCreateSubdomains2D - Creates the index sets for the overlapping Schwarz 15674b9ad928SBarry Smith preconditioner for a two-dimensional problem on a regular grid. 15684b9ad928SBarry Smith 15694b9ad928SBarry Smith Not Collective 15704b9ad928SBarry Smith 15714b9ad928SBarry Smith Input Parameters: 15724b9ad928SBarry Smith + m, n - the number of mesh points in the x and y directions 15734b9ad928SBarry Smith . M, N - the number of subdomains in the x and y directions 15744b9ad928SBarry Smith . dof - degrees of freedom per node 15754b9ad928SBarry Smith - overlap - overlap in mesh lines 15764b9ad928SBarry Smith 15774b9ad928SBarry Smith Output Parameters: 15784b9ad928SBarry Smith + Nsub - the number of subdomains created 15793d03bb48SJed Brown . is - array of index sets defining overlapping (if overlap > 0) subdomains 15803d03bb48SJed Brown - is_local - array of index sets defining non-overlapping subdomains 15814b9ad928SBarry Smith 15824b9ad928SBarry Smith Note: 15834b9ad928SBarry Smith Presently PCAMSCreateSubdomains2d() is valid only for sequential 15844b9ad928SBarry Smith preconditioners. More general related routines are 15854b9ad928SBarry Smith PCASMSetTotalSubdomains() and PCASMSetLocalSubdomains(). 15864b9ad928SBarry Smith 15874b9ad928SBarry Smith Level: advanced 15884b9ad928SBarry Smith 15894b9ad928SBarry Smith .keywords: PC, ASM, additive Schwarz, create, subdomains, 2D, regular grid 15904b9ad928SBarry Smith 15914b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetLocalSubdomains(), PCASMGetSubKSP(), 15924b9ad928SBarry Smith PCASMSetOverlap() 15934b9ad928SBarry Smith @*/ 15947087cfbeSBarry Smith PetscErrorCode PCASMCreateSubdomains2D(PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt dof,PetscInt overlap,PetscInt *Nsub,IS **is,IS **is_local) 15954b9ad928SBarry Smith { 15963d03bb48SJed Brown PetscInt i,j,height,width,ystart,xstart,yleft,yright,xleft,xright,loc_outer; 15976849ba73SBarry Smith PetscErrorCode ierr; 159813f74950SBarry Smith PetscInt nidx,*idx,loc,ii,jj,count; 15994b9ad928SBarry Smith 16004b9ad928SBarry Smith PetscFunctionBegin; 1601e32f2f54SBarry Smith if (dof != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP," "); 16024b9ad928SBarry Smith 16034b9ad928SBarry Smith *Nsub = N*M; 1604854ce69bSBarry Smith ierr = PetscMalloc1(*Nsub,is);CHKERRQ(ierr); 1605854ce69bSBarry Smith ierr = PetscMalloc1(*Nsub,is_local);CHKERRQ(ierr); 16064b9ad928SBarry Smith ystart = 0; 16073d03bb48SJed Brown loc_outer = 0; 16084b9ad928SBarry Smith for (i=0; i<N; i++) { 16094b9ad928SBarry Smith height = n/N + ((n % N) > i); /* height of subdomain */ 1610e32f2f54SBarry Smith if (height < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many N subdomains for mesh dimension n"); 16114b9ad928SBarry Smith yleft = ystart - overlap; if (yleft < 0) yleft = 0; 16124b9ad928SBarry Smith yright = ystart + height + overlap; if (yright > n) yright = n; 16134b9ad928SBarry Smith xstart = 0; 16144b9ad928SBarry Smith for (j=0; j<M; j++) { 16154b9ad928SBarry Smith width = m/M + ((m % M) > j); /* width of subdomain */ 1616e32f2f54SBarry Smith if (width < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many M subdomains for mesh dimension m"); 16174b9ad928SBarry Smith xleft = xstart - overlap; if (xleft < 0) xleft = 0; 16184b9ad928SBarry Smith xright = xstart + width + overlap; if (xright > m) xright = m; 16194b9ad928SBarry Smith nidx = (xright - xleft)*(yright - yleft); 1620785e854fSJed Brown ierr = PetscMalloc1(nidx,&idx);CHKERRQ(ierr); 16214b9ad928SBarry Smith loc = 0; 16224b9ad928SBarry Smith for (ii=yleft; ii<yright; ii++) { 16234b9ad928SBarry Smith count = m*ii + xleft; 16242fa5cd67SKarl Rupp for (jj=xleft; jj<xright; jj++) idx[loc++] = count++; 16254b9ad928SBarry Smith } 162670b3c8c7SBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,nidx,idx,PETSC_COPY_VALUES,(*is)+loc_outer);CHKERRQ(ierr); 16273d03bb48SJed Brown if (overlap == 0) { 16283d03bb48SJed Brown ierr = PetscObjectReference((PetscObject)(*is)[loc_outer]);CHKERRQ(ierr); 16292fa5cd67SKarl Rupp 16303d03bb48SJed Brown (*is_local)[loc_outer] = (*is)[loc_outer]; 16313d03bb48SJed Brown } else { 16323d03bb48SJed Brown for (loc=0,ii=ystart; ii<ystart+height; ii++) { 16333d03bb48SJed Brown for (jj=xstart; jj<xstart+width; jj++) { 16343d03bb48SJed Brown idx[loc++] = m*ii + jj; 16353d03bb48SJed Brown } 16363d03bb48SJed Brown } 163770b3c8c7SBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,loc,idx,PETSC_COPY_VALUES,*is_local+loc_outer);CHKERRQ(ierr); 16383d03bb48SJed Brown } 16394b9ad928SBarry Smith ierr = PetscFree(idx);CHKERRQ(ierr); 16404b9ad928SBarry Smith xstart += width; 16413d03bb48SJed Brown loc_outer++; 16424b9ad928SBarry Smith } 16434b9ad928SBarry Smith ystart += height; 16444b9ad928SBarry Smith } 16454b9ad928SBarry Smith for (i=0; i<*Nsub; i++) { ierr = ISSort((*is)[i]);CHKERRQ(ierr); } 16464b9ad928SBarry Smith PetscFunctionReturn(0); 16474b9ad928SBarry Smith } 16484b9ad928SBarry Smith 16494b9ad928SBarry Smith /*@C 16504b9ad928SBarry Smith PCASMGetLocalSubdomains - Gets the local subdomains (for this processor 16514b9ad928SBarry Smith only) for the additive Schwarz preconditioner. 16524b9ad928SBarry Smith 1653ad4df100SBarry Smith Not Collective 16544b9ad928SBarry Smith 16554b9ad928SBarry Smith Input Parameter: 16564b9ad928SBarry Smith . pc - the preconditioner context 16574b9ad928SBarry Smith 16584b9ad928SBarry Smith Output Parameters: 16594b9ad928SBarry Smith + n - the number of subdomains for this processor (default value = 1) 16602b691e39SMatthew Knepley . is - the index sets that define the subdomains for this processor 16610298fd71SBarry Smith - is_local - the index sets that define the local part of the subdomains for this processor (can be NULL) 16624b9ad928SBarry Smith 16634b9ad928SBarry Smith 16644b9ad928SBarry Smith Notes: 16654b9ad928SBarry Smith The IS numbering is in the parallel, global numbering of the vector. 16664b9ad928SBarry Smith 16674b9ad928SBarry Smith Level: advanced 16684b9ad928SBarry Smith 16694b9ad928SBarry Smith .keywords: PC, ASM, set, local, subdomains, additive Schwarz 16704b9ad928SBarry Smith 16714b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(), 16724b9ad928SBarry Smith PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubmatrices() 16734b9ad928SBarry Smith @*/ 16747087cfbeSBarry Smith PetscErrorCode PCASMGetLocalSubdomains(PC pc,PetscInt *n,IS *is[],IS *is_local[]) 16754b9ad928SBarry Smith { 16762a808120SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 167792bb6962SLisandro Dalcin PetscErrorCode ierr; 1678ace3abfcSBarry Smith PetscBool match; 16794b9ad928SBarry Smith 16804b9ad928SBarry Smith PetscFunctionBegin; 16810700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 16824482741eSBarry Smith PetscValidIntPointer(n,2); 168392bb6962SLisandro Dalcin if (is) PetscValidPointer(is,3); 1684251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)pc,PCASM,&match);CHKERRQ(ierr); 16852a808120SBarry Smith if (!match) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"PC is not a PCASM"); 16864b9ad928SBarry Smith if (n) *n = osm->n_local_true; 16874b9ad928SBarry Smith if (is) *is = osm->is; 16882b691e39SMatthew Knepley if (is_local) *is_local = osm->is_local; 16894b9ad928SBarry Smith PetscFunctionReturn(0); 16904b9ad928SBarry Smith } 16914b9ad928SBarry Smith 16924b9ad928SBarry Smith /*@C 16934b9ad928SBarry Smith PCASMGetLocalSubmatrices - Gets the local submatrices (for this processor 16944b9ad928SBarry Smith only) for the additive Schwarz preconditioner. 16954b9ad928SBarry Smith 1696ad4df100SBarry Smith Not Collective 16974b9ad928SBarry Smith 16984b9ad928SBarry Smith Input Parameter: 16994b9ad928SBarry Smith . pc - the preconditioner context 17004b9ad928SBarry Smith 17014b9ad928SBarry Smith Output Parameters: 17024b9ad928SBarry Smith + n - the number of matrices for this processor (default value = 1) 17034b9ad928SBarry Smith - mat - the matrices 17044b9ad928SBarry Smith 17054b9ad928SBarry Smith Level: advanced 17064b9ad928SBarry Smith 170795452b02SPatrick Sanan Notes: 170895452b02SPatrick Sanan Call after PCSetUp() (or KSPSetUp()) but before PCApply() (or KSPApply()) and before PCSetUpOnBlocks()) 1709cf739d55SBarry Smith 1710cf739d55SBarry Smith Usually one would use PCSetModifySubmatrices() to change the submatrices in building the preconditioner. 1711cf739d55SBarry Smith 17124b9ad928SBarry Smith .keywords: PC, ASM, set, local, subdomains, additive Schwarz, block Jacobi 17134b9ad928SBarry Smith 17144b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(), 1715cf739d55SBarry Smith PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubdomains(), PCSetModifySubmatrices() 17164b9ad928SBarry Smith @*/ 17177087cfbeSBarry Smith PetscErrorCode PCASMGetLocalSubmatrices(PC pc,PetscInt *n,Mat *mat[]) 17184b9ad928SBarry Smith { 17194b9ad928SBarry Smith PC_ASM *osm; 172092bb6962SLisandro Dalcin PetscErrorCode ierr; 1721ace3abfcSBarry Smith PetscBool match; 17224b9ad928SBarry Smith 17234b9ad928SBarry Smith PetscFunctionBegin; 17240700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 172592bb6962SLisandro Dalcin PetscValidIntPointer(n,2); 172692bb6962SLisandro Dalcin if (mat) PetscValidPointer(mat,3); 1727ce94432eSBarry Smith if (!pc->setupcalled) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call after KSPSetUP() or PCSetUp()."); 1728251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)pc,PCASM,&match);CHKERRQ(ierr); 172992bb6962SLisandro Dalcin if (!match) { 173092bb6962SLisandro Dalcin if (n) *n = 0; 17310298fd71SBarry Smith if (mat) *mat = NULL; 173292bb6962SLisandro Dalcin } else { 17334b9ad928SBarry Smith osm = (PC_ASM*)pc->data; 17344b9ad928SBarry Smith if (n) *n = osm->n_local_true; 17354b9ad928SBarry Smith if (mat) *mat = osm->pmat; 173692bb6962SLisandro Dalcin } 17374b9ad928SBarry Smith PetscFunctionReturn(0); 17384b9ad928SBarry Smith } 1739d709ea83SDmitry Karpeev 1740d709ea83SDmitry Karpeev /*@ 1741d709ea83SDmitry Karpeev PCASMSetDMSubdomains - Indicates whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible. 1742f1ee410cSBarry Smith 1743d709ea83SDmitry Karpeev Logically Collective 1744d709ea83SDmitry Karpeev 1745d709ea83SDmitry Karpeev Input Parameter: 1746d709ea83SDmitry Karpeev + pc - the preconditioner 1747d709ea83SDmitry Karpeev - flg - boolean indicating whether to use subdomains defined by the DM 1748d709ea83SDmitry Karpeev 1749d709ea83SDmitry Karpeev Options Database Key: 1750d709ea83SDmitry Karpeev . -pc_asm_dm_subdomains 1751d709ea83SDmitry Karpeev 1752d709ea83SDmitry Karpeev Level: intermediate 1753d709ea83SDmitry Karpeev 1754d709ea83SDmitry Karpeev Notes: 1755d709ea83SDmitry Karpeev PCASMSetTotalSubdomains() and PCASMSetOverlap() take precedence over PCASMSetDMSubdomains(), 1756d709ea83SDmitry Karpeev so setting either of the first two effectively turns the latter off. 1757d709ea83SDmitry Karpeev 1758d709ea83SDmitry Karpeev .keywords: PC, ASM, DM, set, subdomains, additive Schwarz 1759d709ea83SDmitry Karpeev 1760d709ea83SDmitry Karpeev .seealso: PCASMGetDMSubdomains(), PCASMSetTotalSubdomains(), PCASMSetOverlap() 1761d709ea83SDmitry Karpeev PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubdomains() 1762d709ea83SDmitry Karpeev @*/ 1763d709ea83SDmitry Karpeev PetscErrorCode PCASMSetDMSubdomains(PC pc,PetscBool flg) 1764d709ea83SDmitry Karpeev { 1765d709ea83SDmitry Karpeev PC_ASM *osm = (PC_ASM*)pc->data; 1766d709ea83SDmitry Karpeev PetscErrorCode ierr; 1767d709ea83SDmitry Karpeev PetscBool match; 1768d709ea83SDmitry Karpeev 1769d709ea83SDmitry Karpeev PetscFunctionBegin; 1770d709ea83SDmitry Karpeev PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1771d709ea83SDmitry Karpeev PetscValidLogicalCollectiveBool(pc,flg,2); 1772d709ea83SDmitry Karpeev if (pc->setupcalled) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Not for a setup PC."); 1773d709ea83SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)pc,PCASM,&match);CHKERRQ(ierr); 1774d709ea83SDmitry Karpeev if (match) { 1775d709ea83SDmitry Karpeev osm->dm_subdomains = flg; 1776d709ea83SDmitry Karpeev } 1777d709ea83SDmitry Karpeev PetscFunctionReturn(0); 1778d709ea83SDmitry Karpeev } 1779d709ea83SDmitry Karpeev 1780d709ea83SDmitry Karpeev /*@ 1781d709ea83SDmitry Karpeev PCASMGetDMSubdomains - Returns flag indicating whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible. 1782d709ea83SDmitry Karpeev Not Collective 1783d709ea83SDmitry Karpeev 1784d709ea83SDmitry Karpeev Input Parameter: 1785d709ea83SDmitry Karpeev . pc - the preconditioner 1786d709ea83SDmitry Karpeev 1787d709ea83SDmitry Karpeev Output Parameter: 1788d709ea83SDmitry Karpeev . flg - boolean indicating whether to use subdomains defined by the DM 1789d709ea83SDmitry Karpeev 1790d709ea83SDmitry Karpeev Level: intermediate 1791d709ea83SDmitry Karpeev 1792d709ea83SDmitry Karpeev .keywords: PC, ASM, DM, set, subdomains, additive Schwarz 1793d709ea83SDmitry Karpeev 1794d709ea83SDmitry Karpeev .seealso: PCASMSetDMSubdomains(), PCASMSetTotalSubdomains(), PCASMSetOverlap() 1795d709ea83SDmitry Karpeev PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubdomains() 1796d709ea83SDmitry Karpeev @*/ 1797d709ea83SDmitry Karpeev PetscErrorCode PCASMGetDMSubdomains(PC pc,PetscBool* flg) 1798d709ea83SDmitry Karpeev { 1799d709ea83SDmitry Karpeev PC_ASM *osm = (PC_ASM*)pc->data; 1800d709ea83SDmitry Karpeev PetscErrorCode ierr; 1801d709ea83SDmitry Karpeev PetscBool match; 1802d709ea83SDmitry Karpeev 1803d709ea83SDmitry Karpeev PetscFunctionBegin; 1804d709ea83SDmitry Karpeev PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1805d709ea83SDmitry Karpeev PetscValidPointer(flg,2); 1806d709ea83SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)pc,PCASM,&match);CHKERRQ(ierr); 1807d709ea83SDmitry Karpeev if (match) { 1808d709ea83SDmitry Karpeev if (flg) *flg = osm->dm_subdomains; 1809d709ea83SDmitry Karpeev } 1810d709ea83SDmitry Karpeev PetscFunctionReturn(0); 1811d709ea83SDmitry Karpeev } 181280ec0b7dSPatrick Sanan 181380ec0b7dSPatrick Sanan /*@ 181480ec0b7dSPatrick Sanan PCASMGetSubMatType - Gets the matrix type used for ASM subsolves, as a string. 181580ec0b7dSPatrick Sanan 181680ec0b7dSPatrick Sanan Not Collective 181780ec0b7dSPatrick Sanan 181880ec0b7dSPatrick Sanan Input Parameter: 181980ec0b7dSPatrick Sanan . pc - the PC 182080ec0b7dSPatrick Sanan 182180ec0b7dSPatrick Sanan Output Parameter: 1822f1ee410cSBarry Smith . -pc_asm_sub_mat_type - name of matrix type 182380ec0b7dSPatrick Sanan 182480ec0b7dSPatrick Sanan Level: advanced 182580ec0b7dSPatrick Sanan 182680ec0b7dSPatrick Sanan .keywords: PC, PCASM, MatType, set 182780ec0b7dSPatrick Sanan 182880ec0b7dSPatrick Sanan .seealso: PCASMSetSubMatType(), PCASM, PCSetType(), VecSetType(), MatType, Mat 182980ec0b7dSPatrick Sanan @*/ 183080ec0b7dSPatrick Sanan PetscErrorCode PCASMGetSubMatType(PC pc,MatType *sub_mat_type){ 183180ec0b7dSPatrick Sanan PetscErrorCode ierr; 183280ec0b7dSPatrick Sanan 183380ec0b7dSPatrick Sanan ierr = PetscTryMethod(pc,"PCASMGetSubMatType_C",(PC,MatType*),(pc,sub_mat_type));CHKERRQ(ierr); 183480ec0b7dSPatrick Sanan PetscFunctionReturn(0); 183580ec0b7dSPatrick Sanan } 183680ec0b7dSPatrick Sanan 183780ec0b7dSPatrick Sanan /*@ 183880ec0b7dSPatrick Sanan PCASMSetSubMatType - Set the type of matrix used for ASM subsolves 183980ec0b7dSPatrick Sanan 184080ec0b7dSPatrick Sanan Collective on Mat 184180ec0b7dSPatrick Sanan 184280ec0b7dSPatrick Sanan Input Parameters: 184380ec0b7dSPatrick Sanan + pc - the PC object 184480ec0b7dSPatrick Sanan - sub_mat_type - matrix type 184580ec0b7dSPatrick Sanan 184680ec0b7dSPatrick Sanan Options Database Key: 184780ec0b7dSPatrick 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. 184880ec0b7dSPatrick Sanan 184980ec0b7dSPatrick Sanan Notes: 185080ec0b7dSPatrick Sanan See "${PETSC_DIR}/include/petscmat.h" for available types 185180ec0b7dSPatrick Sanan 185280ec0b7dSPatrick Sanan Level: advanced 185380ec0b7dSPatrick Sanan 185480ec0b7dSPatrick Sanan .keywords: PC, PCASM, MatType, set 185580ec0b7dSPatrick Sanan 185680ec0b7dSPatrick Sanan .seealso: PCASMGetSubMatType(), PCASM, PCSetType(), VecSetType(), MatType, Mat 185780ec0b7dSPatrick Sanan @*/ 185880ec0b7dSPatrick Sanan PetscErrorCode PCASMSetSubMatType(PC pc,MatType sub_mat_type) 185980ec0b7dSPatrick Sanan { 186080ec0b7dSPatrick Sanan PetscErrorCode ierr; 186180ec0b7dSPatrick Sanan 186280ec0b7dSPatrick Sanan ierr = PetscTryMethod(pc,"PCASMSetSubMatType_C",(PC,MatType),(pc,sub_mat_type));CHKERRQ(ierr); 186380ec0b7dSPatrick Sanan PetscFunctionReturn(0); 186480ec0b7dSPatrick Sanan } 1865