1dba47a55SKris Buschelman 24b9ad928SBarry Smith /* 34b9ad928SBarry Smith This file defines an additive Schwarz preconditioner for any Mat implementation. 44b9ad928SBarry Smith 54b9ad928SBarry Smith Note that each processor may have any number of subdomains. But in order to 64b9ad928SBarry Smith deal easily with the VecScatter(), we treat each processor as if it has the 74b9ad928SBarry Smith same number of subdomains. 84b9ad928SBarry Smith 94b9ad928SBarry Smith n - total number of true subdomains on all processors 104b9ad928SBarry Smith n_local_true - actual number of subdomains on this processor 114b9ad928SBarry Smith n_local = maximum over all processors of n_local_true 124b9ad928SBarry Smith */ 13af0996ceSBarry Smith #include <petsc/private/pcimpl.h> /*I "petscpc.h" I*/ 141e25c274SJed Brown #include <petscdm.h> 154b9ad928SBarry Smith 164b9ad928SBarry Smith typedef struct { 1713f74950SBarry Smith PetscInt n, n_local, n_local_true; 1813f74950SBarry Smith PetscInt overlap; /* overlap requested by user */ 194b9ad928SBarry Smith KSP *ksp; /* linear solvers for each block */ 20b0de9f23SBarry Smith VecScatter restriction; /* mapping from global to overlapping (process) subdomain*/ 21b0de9f23SBarry Smith VecScatter *lrestriction; /* mapping from subregion to overlapping (process) subdomain */ 22b0de9f23SBarry Smith VecScatter *lprolongation; /* mapping from non-overlapping subregion to overlapping (process) subdomain; used for restrict additive version of algorithms */ 231dd8081eSeaulisa Vec lx, ly; /* work vectors */ 241dd8081eSeaulisa Vec *x,*y; /* work vectors */ 251dd8081eSeaulisa IS lis; /* index set that defines each overlapping multiplicative (process) subdomain */ 262b691e39SMatthew Knepley IS *is; /* index set that defines each overlapping subdomain */ 272b691e39SMatthew Knepley IS *is_local; /* index set that defines each non-overlapping subdomain, may be NULL */ 284b9ad928SBarry Smith Mat *mat,*pmat; /* mat is not currently used */ 294b9ad928SBarry Smith PCASMType type; /* use reduced interpolation, restriction or both */ 30ace3abfcSBarry Smith PetscBool type_set; /* if user set this value (so won't change it for symmetric problems) */ 31ace3abfcSBarry Smith PetscBool same_local_solves; /* flag indicating whether all local solvers are same */ 32ace3abfcSBarry Smith PetscBool sort_indices; /* flag to sort subdomain indices */ 33d709ea83SDmitry Karpeev PetscBool dm_subdomains; /* whether DM is allowed to define subdomains */ 3412cd4985SMatthew G. Knepley PCCompositeType loctype; /* the type of composition for local solves */ 3580ec0b7dSPatrick Sanan MatType sub_mat_type; /* the type of Mat used for subdomain solves (can be MATSAME or NULL) */ 36fb745f2cSMatthew G. Knepley /* For multiplicative solve */ 37235cc4ceSMatthew G. Knepley Mat *lmats; /* submatrices for overlapping multiplicative (process) subdomain */ 384b9ad928SBarry Smith } PC_ASM; 394b9ad928SBarry Smith 406849ba73SBarry Smith static PetscErrorCode PCView_ASM(PC pc,PetscViewer viewer) 414b9ad928SBarry Smith { 4292bb6962SLisandro Dalcin PC_ASM *osm = (PC_ASM*)pc->data; 436849ba73SBarry Smith PetscErrorCode ierr; 4413f74950SBarry Smith PetscMPIInt rank; 454d219a6aSLisandro Dalcin PetscInt i,bsz; 46ace3abfcSBarry Smith PetscBool iascii,isstring; 474b9ad928SBarry Smith PetscViewer sviewer; 484b9ad928SBarry Smith 494b9ad928SBarry Smith PetscFunctionBegin; 50251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 51251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr); 5232077d6dSBarry Smith if (iascii) { 533d03bb48SJed Brown char overlaps[256] = "user-defined overlap",blocks[256] = "total subdomain blocks not yet set"; 548caf3d72SBarry Smith if (osm->overlap >= 0) {ierr = PetscSNPrintf(overlaps,sizeof(overlaps),"amount of overlap = %D",osm->overlap);CHKERRQ(ierr);} 558caf3d72SBarry Smith if (osm->n > 0) {ierr = PetscSNPrintf(blocks,sizeof(blocks),"total subdomain blocks = %D",osm->n);CHKERRQ(ierr);} 56efd4aadfSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %s, %s\n",blocks,overlaps);CHKERRQ(ierr); 57efd4aadfSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," restriction/interpolation type - %s\n",PCASMTypes[osm->type]);CHKERRQ(ierr); 58e64f0791SPatrick Sanan if (osm->dm_subdomains) {ierr = PetscViewerASCIIPrintf(viewer," Additive Schwarz: using DM to define subdomains\n");CHKERRQ(ierr);} 5912cd4985SMatthew G. Knepley if (osm->loctype != PC_COMPOSITE_ADDITIVE) {ierr = PetscViewerASCIIPrintf(viewer," Additive Schwarz: local solve composition type - %s\n",PCCompositeTypes[osm->loctype]);CHKERRQ(ierr);} 60ce94432eSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);CHKERRQ(ierr); 6192bb6962SLisandro Dalcin if (osm->same_local_solves) { 624d219a6aSLisandro Dalcin if (osm->ksp) { 634b9ad928SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Local solve is same for all blocks, in the following KSP and PC objects:\n");CHKERRQ(ierr); 64c5e4d11fSDmitry Karpeev ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 654d219a6aSLisandro Dalcin if (!rank) { 664b9ad928SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 6792bb6962SLisandro Dalcin ierr = KSPView(osm->ksp[0],sviewer);CHKERRQ(ierr); 684b9ad928SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 694b9ad928SBarry Smith } 70c5e4d11fSDmitry Karpeev ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 714d219a6aSLisandro Dalcin } 724b9ad928SBarry Smith } else { 73c5e4d11fSDmitry Karpeev ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 744d219a6aSLisandro Dalcin ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] number of local blocks = %D\n",(int)rank,osm->n_local_true);CHKERRQ(ierr); 754d219a6aSLisandro Dalcin ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 764b9ad928SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Local solve info for each block is in the following KSP and PC objects:\n");CHKERRQ(ierr); 774b9ad928SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 784d219a6aSLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr); 79c5e4d11fSDmitry Karpeev ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 80dd2fa690SBarry Smith for (i=0; i<osm->n_local_true; i++) { 814d219a6aSLisandro Dalcin ierr = ISGetLocalSize(osm->is[i],&bsz);CHKERRQ(ierr); 824d219a6aSLisandro Dalcin ierr = PetscViewerASCIISynchronizedPrintf(sviewer,"[%d] local block number %D, size = %D\n",(int)rank,i,bsz);CHKERRQ(ierr); 8392bb6962SLisandro Dalcin ierr = KSPView(osm->ksp[i],sviewer);CHKERRQ(ierr); 844d219a6aSLisandro Dalcin ierr = PetscViewerASCIISynchronizedPrintf(sviewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr); 854b9ad928SBarry Smith } 86c5e4d11fSDmitry Karpeev ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 874b9ad928SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 884b9ad928SBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 89c5e4d11fSDmitry Karpeev ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 904b9ad928SBarry Smith } 914b9ad928SBarry Smith } else if (isstring) { 924d219a6aSLisandro Dalcin ierr = PetscViewerStringSPrintf(viewer," blocks=%D, overlap=%D, type=%s",osm->n,osm->overlap,PCASMTypes[osm->type]);CHKERRQ(ierr); 93c5e4d11fSDmitry Karpeev ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 9492bb6962SLisandro Dalcin if (osm->ksp) {ierr = KSPView(osm->ksp[0],sviewer);CHKERRQ(ierr);} 95c5e4d11fSDmitry Karpeev ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 964b9ad928SBarry Smith } 974b9ad928SBarry Smith PetscFunctionReturn(0); 984b9ad928SBarry Smith } 994b9ad928SBarry Smith 10092bb6962SLisandro Dalcin static PetscErrorCode PCASMPrintSubdomains(PC pc) 10192bb6962SLisandro Dalcin { 10292bb6962SLisandro Dalcin PC_ASM *osm = (PC_ASM*)pc->data; 10392bb6962SLisandro Dalcin const char *prefix; 10492bb6962SLisandro Dalcin char fname[PETSC_MAX_PATH_LEN+1]; 105643535aeSDmitry Karpeev PetscViewer viewer, sviewer; 106643535aeSDmitry Karpeev char *s; 10792bb6962SLisandro Dalcin PetscInt i,j,nidx; 10892bb6962SLisandro Dalcin const PetscInt *idx; 109643535aeSDmitry Karpeev PetscMPIInt rank, size; 11092bb6962SLisandro Dalcin PetscErrorCode ierr; 111846783a0SBarry Smith 11292bb6962SLisandro Dalcin PetscFunctionBegin; 113ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size);CHKERRQ(ierr); 114ce94432eSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank);CHKERRQ(ierr); 11592bb6962SLisandro Dalcin ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 116c5929fdfSBarry Smith ierr = PetscOptionsGetString(NULL,prefix,"-pc_asm_print_subdomains",fname,PETSC_MAX_PATH_LEN,NULL);CHKERRQ(ierr); 11792bb6962SLisandro Dalcin if (fname[0] == 0) { ierr = PetscStrcpy(fname,"stdout");CHKERRQ(ierr); }; 118ce94432eSBarry Smith ierr = PetscViewerASCIIOpen(PetscObjectComm((PetscObject)pc),fname,&viewer);CHKERRQ(ierr); 119643535aeSDmitry Karpeev for (i=0; i<osm->n_local; i++) { 120643535aeSDmitry Karpeev if (i < osm->n_local_true) { 12192bb6962SLisandro Dalcin ierr = ISGetLocalSize(osm->is[i],&nidx);CHKERRQ(ierr); 12292bb6962SLisandro Dalcin ierr = ISGetIndices(osm->is[i],&idx);CHKERRQ(ierr); 123643535aeSDmitry Karpeev /* Print to a string viewer; no more than 15 characters per index plus 512 char for the header.*/ 124854ce69bSBarry Smith ierr = PetscMalloc1(16*(nidx+1)+512, &s);CHKERRQ(ierr); 125643535aeSDmitry Karpeev ierr = PetscViewerStringOpen(PETSC_COMM_SELF, s, 16*(nidx+1)+512, &sviewer);CHKERRQ(ierr); 126643535aeSDmitry Karpeev ierr = PetscViewerStringSPrintf(sviewer, "[%D:%D] Subdomain %D with overlap:\n", rank, size, i);CHKERRQ(ierr); 12792bb6962SLisandro Dalcin for (j=0; j<nidx; j++) { 128643535aeSDmitry Karpeev ierr = PetscViewerStringSPrintf(sviewer,"%D ",idx[j]);CHKERRQ(ierr); 12992bb6962SLisandro Dalcin } 13092bb6962SLisandro Dalcin ierr = ISRestoreIndices(osm->is[i],&idx);CHKERRQ(ierr); 131643535aeSDmitry Karpeev ierr = PetscViewerStringSPrintf(sviewer,"\n");CHKERRQ(ierr); 132643535aeSDmitry Karpeev ierr = PetscViewerDestroy(&sviewer);CHKERRQ(ierr); 133c5e4d11fSDmitry Karpeev ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 134643535aeSDmitry Karpeev ierr = PetscViewerASCIISynchronizedPrintf(viewer, s);CHKERRQ(ierr); 135643535aeSDmitry Karpeev ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 136c5e4d11fSDmitry Karpeev ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 137643535aeSDmitry Karpeev ierr = PetscFree(s);CHKERRQ(ierr); 1382b691e39SMatthew Knepley if (osm->is_local) { 139643535aeSDmitry Karpeev /* Print to a string viewer; no more than 15 characters per index plus 512 char for the header.*/ 140854ce69bSBarry Smith ierr = PetscMalloc1(16*(nidx+1)+512, &s);CHKERRQ(ierr); 141643535aeSDmitry Karpeev ierr = PetscViewerStringOpen(PETSC_COMM_SELF, s, 16*(nidx+1)+512, &sviewer);CHKERRQ(ierr); 14209d011a0SDmitry Karpeev ierr = PetscViewerStringSPrintf(sviewer, "[%D:%D] Subdomain %D without overlap:\n", rank, size, i);CHKERRQ(ierr); 1432b691e39SMatthew Knepley ierr = ISGetLocalSize(osm->is_local[i],&nidx);CHKERRQ(ierr); 1442b691e39SMatthew Knepley ierr = ISGetIndices(osm->is_local[i],&idx);CHKERRQ(ierr); 1452b691e39SMatthew Knepley for (j=0; j<nidx; j++) { 14609d011a0SDmitry Karpeev ierr = PetscViewerStringSPrintf(sviewer,"%D ",idx[j]);CHKERRQ(ierr); 1472b691e39SMatthew Knepley } 1482b691e39SMatthew Knepley ierr = ISRestoreIndices(osm->is_local[i],&idx);CHKERRQ(ierr); 14909d011a0SDmitry Karpeev ierr = PetscViewerStringSPrintf(sviewer,"\n");CHKERRQ(ierr); 150643535aeSDmitry Karpeev ierr = PetscViewerDestroy(&sviewer);CHKERRQ(ierr); 151c5e4d11fSDmitry Karpeev ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 152643535aeSDmitry Karpeev ierr = PetscViewerASCIISynchronizedPrintf(viewer, s);CHKERRQ(ierr); 153643535aeSDmitry Karpeev ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 154c5e4d11fSDmitry Karpeev ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 155643535aeSDmitry Karpeev ierr = PetscFree(s);CHKERRQ(ierr); 156643535aeSDmitry Karpeev } 1572fa5cd67SKarl Rupp } else { 158643535aeSDmitry Karpeev /* Participate in collective viewer calls. */ 159c5e4d11fSDmitry Karpeev ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 160643535aeSDmitry Karpeev ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 161c5e4d11fSDmitry Karpeev ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 162643535aeSDmitry Karpeev /* Assume either all ranks have is_local or none do. */ 163643535aeSDmitry Karpeev if (osm->is_local) { 164c5e4d11fSDmitry Karpeev ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 165643535aeSDmitry Karpeev ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 166c5e4d11fSDmitry Karpeev ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 167643535aeSDmitry Karpeev } 168fdc48646SDmitry Karpeev } 16992bb6962SLisandro Dalcin } 17092bb6962SLisandro Dalcin ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 171fcfd50ebSBarry Smith ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 17292bb6962SLisandro Dalcin PetscFunctionReturn(0); 17392bb6962SLisandro Dalcin } 17492bb6962SLisandro Dalcin 1756849ba73SBarry Smith static PetscErrorCode PCSetUp_ASM(PC pc) 1764b9ad928SBarry Smith { 1774b9ad928SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 1786849ba73SBarry Smith PetscErrorCode ierr; 179ace3abfcSBarry Smith PetscBool symset,flg; 18087e86712SBarry Smith PetscInt i,m,m_local; 1814b9ad928SBarry Smith MatReuse scall = MAT_REUSE_MATRIX; 1824b9ad928SBarry Smith IS isl; 1834b9ad928SBarry Smith KSP ksp; 1844b9ad928SBarry Smith PC subpc; 1852dcb1b2aSMatthew Knepley const char *prefix,*pprefix; 18623ce1328SBarry Smith Vec vec; 1870298fd71SBarry Smith DM *domain_dm = NULL; 1884b9ad928SBarry Smith 1894b9ad928SBarry Smith PetscFunctionBegin; 1904b9ad928SBarry Smith if (!pc->setupcalled) { 191265a2120SBarry Smith PetscInt m; 19292bb6962SLisandro Dalcin 19392bb6962SLisandro Dalcin if (!osm->type_set) { 19492bb6962SLisandro Dalcin ierr = MatIsSymmetricKnown(pc->pmat,&symset,&flg);CHKERRQ(ierr); 1952fa5cd67SKarl Rupp if (symset && flg) osm->type = PC_ASM_BASIC; 19692bb6962SLisandro Dalcin } 19792bb6962SLisandro Dalcin 1982b837212SDmitry Karpeev /* Note: if subdomains have been set either via PCASMSetTotalSubdomains() or via PCASMSetLocalSubdomains(), osm->n_local_true will not be PETSC_DECIDE */ 1992b837212SDmitry Karpeev if (osm->n_local_true == PETSC_DECIDE) { 20069ca1f37SDmitry Karpeev /* no subdomains given */ 20165db9045SDmitry Karpeev /* try pc->dm first, if allowed */ 202d709ea83SDmitry Karpeev if (osm->dm_subdomains && pc->dm) { 203feb221f8SDmitry Karpeev PetscInt num_domains, d; 204feb221f8SDmitry Karpeev char **domain_names; 2058d4ac253SDmitry Karpeev IS *inner_domain_is, *outer_domain_is; 2068d4ac253SDmitry Karpeev ierr = DMCreateDomainDecomposition(pc->dm, &num_domains, &domain_names, &inner_domain_is, &outer_domain_is, &domain_dm);CHKERRQ(ierr); 207704f0395SPatrick Sanan osm->overlap = -1; /* We do not want to increase the overlap of the IS. 208704f0395SPatrick Sanan A future improvement of this code might allow one to use 209704f0395SPatrick Sanan DM-defined subdomains and also increase the overlap, 210704f0395SPatrick Sanan but that is not currently supported */ 21169ca1f37SDmitry Karpeev if (num_domains) { 2128d4ac253SDmitry Karpeev ierr = PCASMSetLocalSubdomains(pc, num_domains, outer_domain_is, inner_domain_is);CHKERRQ(ierr); 21369ca1f37SDmitry Karpeev } 214feb221f8SDmitry Karpeev for (d = 0; d < num_domains; ++d) { 215a35b7b57SDmitry Karpeev if (domain_names) {ierr = PetscFree(domain_names[d]);CHKERRQ(ierr);} 216a35b7b57SDmitry Karpeev if (inner_domain_is) {ierr = ISDestroy(&inner_domain_is[d]);CHKERRQ(ierr);} 217a35b7b57SDmitry Karpeev if (outer_domain_is) {ierr = ISDestroy(&outer_domain_is[d]);CHKERRQ(ierr);} 218feb221f8SDmitry Karpeev } 219feb221f8SDmitry Karpeev ierr = PetscFree(domain_names);CHKERRQ(ierr); 2208d4ac253SDmitry Karpeev ierr = PetscFree(inner_domain_is);CHKERRQ(ierr); 2218d4ac253SDmitry Karpeev ierr = PetscFree(outer_domain_is);CHKERRQ(ierr); 222feb221f8SDmitry Karpeev } 2232b837212SDmitry Karpeev if (osm->n_local_true == PETSC_DECIDE) { 22469ca1f37SDmitry Karpeev /* still no subdomains; use one subdomain per processor */ 2252b837212SDmitry Karpeev osm->n_local_true = 1; 226feb221f8SDmitry Karpeev } 2272b837212SDmitry Karpeev } 2282b837212SDmitry Karpeev { /* determine the global and max number of subdomains */ 2296ac3741eSJed Brown struct {PetscInt max,sum;} inwork,outwork; 230c10200c1SHong Zhang PetscMPIInt size; 231c10200c1SHong Zhang 2326ac3741eSJed Brown inwork.max = osm->n_local_true; 2336ac3741eSJed Brown inwork.sum = osm->n_local_true; 234367daffbSBarry Smith ierr = MPIU_Allreduce(&inwork,&outwork,1,MPIU_2INT,MPIU_MAXSUM_OP,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr); 2356ac3741eSJed Brown osm->n_local = outwork.max; 2366ac3741eSJed Brown osm->n = outwork.sum; 237c10200c1SHong Zhang 238c10200c1SHong Zhang ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr); 239c10200c1SHong Zhang if (outwork.max == 1 && outwork.sum == size) { 2407dae84e0SHong Zhang /* osm->n_local_true = 1 on all processes, set this option may enable use of optimized MatCreateSubMatrices() implementation */ 241c10200c1SHong Zhang ierr = MatSetOption(pc->pmat,MAT_SUBMAT_SINGLEIS,PETSC_TRUE);CHKERRQ(ierr); 242c10200c1SHong Zhang } 2434b9ad928SBarry Smith } 24478904715SLisandro Dalcin if (!osm->is) { /* create the index sets */ 24578904715SLisandro Dalcin ierr = PCASMCreateSubdomains(pc->pmat,osm->n_local_true,&osm->is);CHKERRQ(ierr); 2464b9ad928SBarry Smith } 247f5234e35SJed Brown if (osm->n_local_true > 1 && !osm->is_local) { 248785e854fSJed Brown ierr = PetscMalloc1(osm->n_local_true,&osm->is_local);CHKERRQ(ierr); 249f5234e35SJed Brown for (i=0; i<osm->n_local_true; i++) { 250f5234e35SJed Brown if (osm->overlap > 0) { /* With positive overlap, osm->is[i] will be modified */ 251f5234e35SJed Brown ierr = ISDuplicate(osm->is[i],&osm->is_local[i]);CHKERRQ(ierr); 252f5234e35SJed Brown ierr = ISCopy(osm->is[i],osm->is_local[i]);CHKERRQ(ierr); 253f5234e35SJed Brown } else { 254f5234e35SJed Brown ierr = PetscObjectReference((PetscObject)osm->is[i]);CHKERRQ(ierr); 255f5234e35SJed Brown osm->is_local[i] = osm->is[i]; 256f5234e35SJed Brown } 257f5234e35SJed Brown } 258f5234e35SJed Brown } 25992bb6962SLisandro Dalcin ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 26090d69ab7SBarry Smith flg = PETSC_FALSE; 261c5929fdfSBarry Smith ierr = PetscOptionsGetBool(NULL,prefix,"-pc_asm_print_subdomains",&flg,NULL);CHKERRQ(ierr); 26292bb6962SLisandro Dalcin if (flg) { ierr = PCASMPrintSubdomains(pc);CHKERRQ(ierr); } 2634b9ad928SBarry Smith 2643d03bb48SJed Brown if (osm->overlap > 0) { 2654b9ad928SBarry Smith /* Extend the "overlapping" regions by a number of steps */ 26692bb6962SLisandro Dalcin ierr = MatIncreaseOverlap(pc->pmat,osm->n_local_true,osm->is,osm->overlap);CHKERRQ(ierr); 2673d03bb48SJed Brown } 2686ed231c7SMatthew Knepley if (osm->sort_indices) { 26992bb6962SLisandro Dalcin for (i=0; i<osm->n_local_true; i++) { 2704b9ad928SBarry Smith ierr = ISSort(osm->is[i]);CHKERRQ(ierr); 2712b691e39SMatthew Knepley if (osm->is_local) { 2722b691e39SMatthew Knepley ierr = ISSort(osm->is_local[i]);CHKERRQ(ierr); 2732b691e39SMatthew Knepley } 2744b9ad928SBarry Smith } 2756ed231c7SMatthew Knepley } 2764b9ad928SBarry Smith 277f6991133SBarry Smith if (!osm->ksp) { 27878904715SLisandro Dalcin /* Create the local solvers */ 279785e854fSJed Brown ierr = PetscMalloc1(osm->n_local_true,&osm->ksp);CHKERRQ(ierr); 280feb221f8SDmitry Karpeev if (domain_dm) { 281feb221f8SDmitry Karpeev ierr = PetscInfo(pc,"Setting up ASM subproblems using the embedded DM\n");CHKERRQ(ierr); 282feb221f8SDmitry Karpeev } 28392bb6962SLisandro Dalcin for (i=0; i<osm->n_local_true; i++) { 2844b9ad928SBarry Smith ierr = KSPCreate(PETSC_COMM_SELF,&ksp);CHKERRQ(ierr); 285422a814eSBarry Smith ierr = KSPSetErrorIfNotConverged(ksp,pc->erroriffailure);CHKERRQ(ierr); 2863bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ksp);CHKERRQ(ierr); 28792bb6962SLisandro Dalcin ierr = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,1);CHKERRQ(ierr); 2884b9ad928SBarry Smith ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr); 2894b9ad928SBarry Smith ierr = KSPGetPC(ksp,&subpc);CHKERRQ(ierr); 2904b9ad928SBarry Smith ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 2914b9ad928SBarry Smith ierr = KSPSetOptionsPrefix(ksp,prefix);CHKERRQ(ierr); 2924b9ad928SBarry Smith ierr = KSPAppendOptionsPrefix(ksp,"sub_");CHKERRQ(ierr); 293feb221f8SDmitry Karpeev if (domain_dm) { 294feb221f8SDmitry Karpeev ierr = KSPSetDM(ksp, domain_dm[i]);CHKERRQ(ierr); 295feb221f8SDmitry Karpeev ierr = KSPSetDMActive(ksp, PETSC_FALSE);CHKERRQ(ierr); 296feb221f8SDmitry Karpeev ierr = DMDestroy(&domain_dm[i]);CHKERRQ(ierr); 297feb221f8SDmitry Karpeev } 2984b9ad928SBarry Smith osm->ksp[i] = ksp; 2994b9ad928SBarry Smith } 300feb221f8SDmitry Karpeev if (domain_dm) { 301feb221f8SDmitry Karpeev ierr = PetscFree(domain_dm);CHKERRQ(ierr); 302feb221f8SDmitry Karpeev } 303f6991133SBarry Smith } 3041dd8081eSeaulisa 305fb745f2cSMatthew G. Knepley ierr = ISConcatenate(PETSC_COMM_SELF, osm->n_local_true, osm->is, &osm->lis);CHKERRQ(ierr); 306fb745f2cSMatthew G. Knepley ierr = ISSortRemoveDups(osm->lis);CHKERRQ(ierr); 307fb745f2cSMatthew G. Knepley ierr = ISGetLocalSize(osm->lis, &m);CHKERRQ(ierr); 308fb745f2cSMatthew G. Knepley ierr = VecCreateSeq(PETSC_COMM_SELF, m, &osm->lx);CHKERRQ(ierr); 309fb745f2cSMatthew G. Knepley ierr = VecDuplicate(osm->lx, &osm->ly);CHKERRQ(ierr); 3101dd8081eSeaulisa 3114b9ad928SBarry Smith scall = MAT_INITIAL_MATRIX; 3124b9ad928SBarry Smith } else { 3134b9ad928SBarry Smith /* 3144b9ad928SBarry Smith Destroy the blocks from the previous iteration 3154b9ad928SBarry Smith */ 316e09e08ccSBarry Smith if (pc->flag == DIFFERENT_NONZERO_PATTERN) { 3174b9ad928SBarry Smith ierr = MatDestroyMatrices(osm->n_local_true,&osm->pmat);CHKERRQ(ierr); 3184b9ad928SBarry Smith scall = MAT_INITIAL_MATRIX; 3194b9ad928SBarry Smith } 3204b9ad928SBarry Smith } 3214b9ad928SBarry Smith 32292bb6962SLisandro Dalcin /* 32392bb6962SLisandro Dalcin Extract out the submatrices 32492bb6962SLisandro Dalcin */ 3257dae84e0SHong Zhang ierr = MatCreateSubMatrices(pc->pmat,osm->n_local_true,osm->is,osm->is,scall,&osm->pmat);CHKERRQ(ierr); 32692bb6962SLisandro Dalcin if (scall == MAT_INITIAL_MATRIX) { 32792bb6962SLisandro Dalcin ierr = PetscObjectGetOptionsPrefix((PetscObject)pc->pmat,&pprefix);CHKERRQ(ierr); 32892bb6962SLisandro Dalcin for (i=0; i<osm->n_local_true; i++) { 3293bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)osm->pmat[i]);CHKERRQ(ierr); 33078904715SLisandro Dalcin ierr = PetscObjectSetOptionsPrefix((PetscObject)osm->pmat[i],pprefix);CHKERRQ(ierr); 33192bb6962SLisandro Dalcin } 33292bb6962SLisandro Dalcin } 33380ec0b7dSPatrick Sanan 33480ec0b7dSPatrick Sanan /* Convert the types of the submatrices (if needbe) */ 33580ec0b7dSPatrick Sanan if (osm->sub_mat_type) { 33680ec0b7dSPatrick Sanan for (i=0; i<osm->n_local_true; i++) { 33780ec0b7dSPatrick Sanan ierr = MatConvert(osm->pmat[i],osm->sub_mat_type,MAT_INPLACE_MATRIX,&(osm->pmat[i]));CHKERRQ(ierr); 33880ec0b7dSPatrick Sanan } 33980ec0b7dSPatrick Sanan } 34080ec0b7dSPatrick Sanan 34180ec0b7dSPatrick Sanan if(!pc->setupcalled){ 34280ec0b7dSPatrick Sanan /* Create the local work vectors (from the local matrices) and scatter contexts */ 34380ec0b7dSPatrick Sanan ierr = MatCreateVecs(pc->pmat,&vec,0);CHKERRQ(ierr); 3445b3c0d42Seaulisa 3451dd8081eSeaulisa 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()"); 3461dd8081eSeaulisa if (osm->is_local && osm->type == PC_ASM_RESTRICT && osm->loctype == PC_COMPOSITE_ADDITIVE) { 3471dd8081eSeaulisa ierr = PetscMalloc1(osm->n_local_true,&osm->lprolongation);CHKERRQ(ierr); 3481dd8081eSeaulisa } 3491dd8081eSeaulisa ierr = PetscMalloc1(osm->n_local_true,&osm->lrestriction);CHKERRQ(ierr); 3501dd8081eSeaulisa ierr = PetscMalloc1(osm->n_local_true,&osm->x);CHKERRQ(ierr); 3511dd8081eSeaulisa ierr = PetscMalloc1(osm->n_local_true,&osm->y);CHKERRQ(ierr); 352347574c9Seaulisa 353347574c9Seaulisa ierr = ISGetLocalSize(osm->lis,&m);CHKERRQ(ierr); 354347574c9Seaulisa ierr = ISCreateStride(PETSC_COMM_SELF,m,0,1,&isl);CHKERRQ(ierr); 3551dd8081eSeaulisa ierr = VecScatterCreate(vec,osm->lis,osm->lx,isl,&osm->restriction);CHKERRQ(ierr); 356347574c9Seaulisa ierr = ISDestroy(&isl);CHKERRQ(ierr); 357347574c9Seaulisa 358347574c9Seaulisa 35980ec0b7dSPatrick Sanan for (i=0; i<osm->n_local_true; ++i) { 3605b3c0d42Seaulisa ISLocalToGlobalMapping ltog; 3615b3c0d42Seaulisa IS isll; 3625b3c0d42Seaulisa const PetscInt *idx_is; 3635b3c0d42Seaulisa PetscInt *idx_lis,nout; 3645b3c0d42Seaulisa 36555817e79SBarry Smith ierr = ISGetLocalSize(osm->is[i],&m);CHKERRQ(ierr); 36655817e79SBarry Smith ierr = MatCreateVecs(osm->pmat[i],&osm->x[i],NULL);CHKERRQ(ierr); 36755817e79SBarry Smith ierr = VecDuplicate(osm->x[i],&osm->y[i]);CHKERRQ(ierr); 36855817e79SBarry Smith 369b0de9f23SBarry Smith /* generate a scatter from ly to y[i] picking all the overlapping is[i] entries */ 3705b3c0d42Seaulisa ierr = ISLocalToGlobalMappingCreateIS(osm->lis,<og);CHKERRQ(ierr); 3715b3c0d42Seaulisa ierr = ISGetLocalSize(osm->is[i],&m);CHKERRQ(ierr); 3725b3c0d42Seaulisa ierr = ISGetIndices(osm->is[i], &idx_is);CHKERRQ(ierr); 3735b3c0d42Seaulisa ierr = PetscMalloc1(m,&idx_lis);CHKERRQ(ierr); 3745b3c0d42Seaulisa ierr = ISGlobalToLocalMappingApply(ltog,IS_GTOLM_DROP,m,idx_is,&nout,idx_lis);CHKERRQ(ierr); 3755b3c0d42Seaulisa if (nout != m) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"is not a subset of lis"); 3765b3c0d42Seaulisa ierr = ISRestoreIndices(osm->is[i], &idx_is);CHKERRQ(ierr); 3775b3c0d42Seaulisa ierr = ISCreateGeneral(PETSC_COMM_SELF,m,idx_lis,PETSC_OWN_POINTER,&isll);CHKERRQ(ierr); 378d8b3b5e3Seaulisa ierr = ISLocalToGlobalMappingDestroy(<og);CHKERRQ(ierr); 3795b3c0d42Seaulisa ierr = ISCreateStride(PETSC_COMM_SELF,m,0,1,&isl);CHKERRQ(ierr); 3801dd8081eSeaulisa ierr = VecScatterCreate(osm->ly,isll,osm->y[i],isl,&osm->lrestriction[i]);CHKERRQ(ierr); 3815b3c0d42Seaulisa ierr = ISDestroy(&isll);CHKERRQ(ierr); 3825b3c0d42Seaulisa ierr = ISDestroy(&isl);CHKERRQ(ierr); 383b0de9f23SBarry Smith if (osm->lprolongation) { /* generate a scatter from y[i] to ly picking only the the non-overalapping is_local[i] entries */ 384d8b3b5e3Seaulisa ISLocalToGlobalMapping ltog; 385d8b3b5e3Seaulisa IS isll,isll_local; 386d8b3b5e3Seaulisa const PetscInt *idx_local; 387d8b3b5e3Seaulisa PetscInt *idx1, *idx2, nout; 388d8b3b5e3Seaulisa 389d8b3b5e3Seaulisa ierr = ISGetLocalSize(osm->is_local[i],&m_local);CHKERRQ(ierr); 390d8b3b5e3Seaulisa ierr = ISGetIndices(osm->is_local[i], &idx_local);CHKERRQ(ierr); 391d8b3b5e3Seaulisa 392d8b3b5e3Seaulisa ierr = ISLocalToGlobalMappingCreateIS(osm->is[i],<og);CHKERRQ(ierr); 393d8b3b5e3Seaulisa ierr = PetscMalloc1(m_local,&idx1);CHKERRQ(ierr); 394d8b3b5e3Seaulisa ierr = ISGlobalToLocalMappingApply(ltog,IS_GTOLM_DROP,m_local,idx_local,&nout,idx1);CHKERRQ(ierr); 395d8b3b5e3Seaulisa ierr = ISLocalToGlobalMappingDestroy(<og);CHKERRQ(ierr); 396d8b3b5e3Seaulisa if (nout != m_local) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"is_local not a subset of is"); 397d8b3b5e3Seaulisa ierr = ISCreateGeneral(PETSC_COMM_SELF,m_local,idx1,PETSC_OWN_POINTER,&isll);CHKERRQ(ierr); 398d8b3b5e3Seaulisa 399d8b3b5e3Seaulisa ierr = ISLocalToGlobalMappingCreateIS(osm->lis,<og);CHKERRQ(ierr); 400d8b3b5e3Seaulisa ierr = PetscMalloc1(m_local,&idx2);CHKERRQ(ierr); 401d8b3b5e3Seaulisa ierr = ISGlobalToLocalMappingApply(ltog,IS_GTOLM_DROP,m_local,idx_local,&nout,idx2);CHKERRQ(ierr); 402d8b3b5e3Seaulisa ierr = ISLocalToGlobalMappingDestroy(<og);CHKERRQ(ierr); 403d8b3b5e3Seaulisa if (nout != m_local) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"is_local not a subset of lis"); 404d8b3b5e3Seaulisa ierr = ISCreateGeneral(PETSC_COMM_SELF,m_local,idx2,PETSC_OWN_POINTER,&isll_local);CHKERRQ(ierr); 405d8b3b5e3Seaulisa 406d8b3b5e3Seaulisa ierr = ISRestoreIndices(osm->is_local[i], &idx_local);CHKERRQ(ierr); 4071dd8081eSeaulisa ierr = VecScatterCreate(osm->y[i],isll,osm->ly,isll_local,&osm->lprolongation[i]);CHKERRQ(ierr); 408d8b3b5e3Seaulisa 409d8b3b5e3Seaulisa ierr = ISDestroy(&isll);CHKERRQ(ierr); 410d8b3b5e3Seaulisa ierr = ISDestroy(&isll_local);CHKERRQ(ierr); 411d8b3b5e3Seaulisa } 41280ec0b7dSPatrick Sanan } 41380ec0b7dSPatrick Sanan ierr = VecDestroy(&vec);CHKERRQ(ierr); 41480ec0b7dSPatrick Sanan } 41580ec0b7dSPatrick Sanan 416fb745f2cSMatthew G. Knepley if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE) { 417235cc4ceSMatthew G. Knepley IS *cis; 418235cc4ceSMatthew G. Knepley PetscInt c; 419235cc4ceSMatthew G. Knepley 420235cc4ceSMatthew G. Knepley ierr = PetscMalloc1(osm->n_local_true, &cis);CHKERRQ(ierr); 421235cc4ceSMatthew G. Knepley for (c = 0; c < osm->n_local_true; ++c) cis[c] = osm->lis; 4227dae84e0SHong Zhang ierr = MatCreateSubMatrices(pc->pmat, osm->n_local_true, osm->is, cis, scall, &osm->lmats);CHKERRQ(ierr); 423235cc4ceSMatthew G. Knepley ierr = PetscFree(cis);CHKERRQ(ierr); 424fb745f2cSMatthew G. Knepley } 4254b9ad928SBarry Smith 4264b9ad928SBarry Smith /* Return control to the user so that the submatrices can be modified (e.g., to apply 4274b9ad928SBarry Smith different boundary conditions for the submatrices than for the global problem) */ 42892bb6962SLisandro Dalcin ierr = PCModifySubMatrices(pc,osm->n_local_true,osm->is,osm->is,osm->pmat,pc->modifysubmatricesP);CHKERRQ(ierr); 4294b9ad928SBarry Smith 43092bb6962SLisandro Dalcin /* 43192bb6962SLisandro Dalcin Loop over subdomains putting them into local ksp 43292bb6962SLisandro Dalcin */ 43392bb6962SLisandro Dalcin for (i=0; i<osm->n_local_true; i++) { 43423ee1639SBarry Smith ierr = KSPSetOperators(osm->ksp[i],osm->pmat[i],osm->pmat[i]);CHKERRQ(ierr); 43592bb6962SLisandro Dalcin if (!pc->setupcalled) { 436bf108f30SBarry Smith ierr = KSPSetFromOptions(osm->ksp[i]);CHKERRQ(ierr); 4374b9ad928SBarry Smith } 43892bb6962SLisandro Dalcin } 4394b9ad928SBarry Smith PetscFunctionReturn(0); 4404b9ad928SBarry Smith } 4414b9ad928SBarry Smith 4426849ba73SBarry Smith static PetscErrorCode PCSetUpOnBlocks_ASM(PC pc) 4434b9ad928SBarry Smith { 4444b9ad928SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 4456849ba73SBarry Smith PetscErrorCode ierr; 44613f74950SBarry Smith PetscInt i; 447539c167fSBarry Smith KSPConvergedReason reason; 4484b9ad928SBarry Smith 4494b9ad928SBarry Smith PetscFunctionBegin; 4504b9ad928SBarry Smith for (i=0; i<osm->n_local_true; i++) { 4514b9ad928SBarry Smith ierr = KSPSetUp(osm->ksp[i]);CHKERRQ(ierr); 452539c167fSBarry Smith ierr = KSPGetConvergedReason(osm->ksp[i],&reason);CHKERRQ(ierr); 453539c167fSBarry Smith if (reason == KSP_DIVERGED_PCSETUP_FAILED) { 454261222b2SHong Zhang pc->failedreason = PC_SUBPC_ERROR; 455e0eafd54SHong Zhang } 4564b9ad928SBarry Smith } 4574b9ad928SBarry Smith PetscFunctionReturn(0); 4584b9ad928SBarry Smith } 4594b9ad928SBarry Smith 4606849ba73SBarry Smith static PetscErrorCode PCApply_ASM(PC pc,Vec x,Vec y) 4614b9ad928SBarry Smith { 4624b9ad928SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 4636849ba73SBarry Smith PetscErrorCode ierr; 4641dd8081eSeaulisa PetscInt i,n_local_true = osm->n_local_true; 4654b9ad928SBarry Smith ScatterMode forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE; 4664b9ad928SBarry Smith 4674b9ad928SBarry Smith PetscFunctionBegin; 4684b9ad928SBarry Smith /* 4694b9ad928SBarry Smith Support for limiting the restriction or interpolation to only local 4704b9ad928SBarry Smith subdomain values (leaving the other values 0). 4714b9ad928SBarry Smith */ 4724b9ad928SBarry Smith if (!(osm->type & PC_ASM_RESTRICT)) { 4734b9ad928SBarry Smith forward = SCATTER_FORWARD_LOCAL; 4744b9ad928SBarry Smith /* have to zero the work RHS since scatter may leave some slots empty */ 4751dd8081eSeaulisa ierr = VecSet(osm->lx, 0.0);CHKERRQ(ierr); 4764b9ad928SBarry Smith } 477347574c9Seaulisa if (!(osm->type & PC_ASM_INTERPOLATE)) { 478347574c9Seaulisa reverse = SCATTER_REVERSE_LOCAL; 479347574c9Seaulisa } 4804b9ad928SBarry Smith 481347574c9Seaulisa if(osm->loctype == PC_COMPOSITE_MULTIPLICATIVE || osm->loctype == PC_COMPOSITE_ADDITIVE){ 482b0de9f23SBarry Smith /* zero the global and the local solutions */ 48312cd4985SMatthew G. Knepley ierr = VecZeroEntries(y);CHKERRQ(ierr); 4845b3c0d42Seaulisa ierr = VecSet(osm->ly, 0.0);CHKERRQ(ierr); 485347574c9Seaulisa 486b0de9f23SBarry Smith /* Copy the global RHS to local RHS including the ghost nodes */ 4871dd8081eSeaulisa ierr = VecScatterBegin(osm->restriction, x, osm->lx, INSERT_VALUES, forward);CHKERRQ(ierr); 4881dd8081eSeaulisa ierr = VecScatterEnd(osm->restriction, x, osm->lx, INSERT_VALUES, forward);CHKERRQ(ierr); 489347574c9Seaulisa 490b0de9f23SBarry Smith /* Restrict local RHS to the overlapping 0-block RHS */ 4911dd8081eSeaulisa ierr = VecScatterBegin(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward);CHKERRQ(ierr); 4921dd8081eSeaulisa ierr = VecScatterEnd(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward);CHKERRQ(ierr); 493d8b3b5e3Seaulisa 49412cd4985SMatthew G. Knepley /* do the local solves */ 49512cd4985SMatthew G. Knepley for (i = 0; i < n_local_true; ++i) { 496347574c9Seaulisa 497b0de9f23SBarry Smith /* solve the overlapping i-block */ 49812cd4985SMatthew G. Knepley ierr = KSPSolve(osm->ksp[i], osm->x[i], osm->y[i]);CHKERRQ(ierr); 499d8b3b5e3Seaulisa 500b0de9f23SBarry Smith if (osm->lprolongation) { /* interpolate the non-overalapping i-block solution to the local solution (only for restrictive additive) */ 5011dd8081eSeaulisa ierr = VecScatterBegin(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward);CHKERRQ(ierr); 5021dd8081eSeaulisa ierr = VecScatterEnd(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward);CHKERRQ(ierr); 503d8b3b5e3Seaulisa } 504b0de9f23SBarry Smith else{ /* interpolate the overalapping i-block solution to the local solution */ 5051dd8081eSeaulisa ierr = VecScatterBegin(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse);CHKERRQ(ierr); 5061dd8081eSeaulisa ierr = VecScatterEnd(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse);CHKERRQ(ierr); 507d8b3b5e3Seaulisa } 508347574c9Seaulisa 509347574c9Seaulisa if (i < n_local_true-1) { 510b0de9f23SBarry Smith /* Restrict local RHS to the overlapping (i+1)-block RHS */ 5111dd8081eSeaulisa ierr = VecScatterBegin(osm->lrestriction[i+1], osm->lx, osm->x[i+1], INSERT_VALUES, forward);CHKERRQ(ierr); 5121dd8081eSeaulisa ierr = VecScatterEnd(osm->lrestriction[i+1], osm->lx, osm->x[i+1], INSERT_VALUES, forward);CHKERRQ(ierr); 513347574c9Seaulisa 514347574c9Seaulisa if ( osm->loctype == PC_COMPOSITE_MULTIPLICATIVE){ 515b0de9f23SBarry Smith /* udpdate the overlapping (i+1)-block RHS using the current local solution */ 516347574c9Seaulisa ierr = MatMult(osm->lmats[i+1], osm->ly, osm->y[i+1]);CHKERRQ(ierr); 517347574c9Seaulisa ierr = VecAXPBY(osm->x[i+1],-1.,1., osm->y[i+1]); CHKERRQ(ierr); 5187c3d802fSMatthew G. Knepley } 51912cd4985SMatthew G. Knepley } 52012cd4985SMatthew G. Knepley } 521b0de9f23SBarry Smith /* Add the local solution to the global solution including the ghost nodes */ 5221dd8081eSeaulisa ierr = VecScatterBegin(osm->restriction, osm->ly, y, ADD_VALUES, reverse);CHKERRQ(ierr); 5231dd8081eSeaulisa ierr = VecScatterEnd(osm->restriction, osm->ly, y, ADD_VALUES, reverse);CHKERRQ(ierr); 524347574c9Seaulisa }else{ 525347574c9Seaulisa SETERRQ1(PetscObjectComm((PetscObject) pc), PETSC_ERR_ARG_WRONG, "Invalid local composition type: %s", PCCompositeTypes[osm->loctype]); 52612cd4985SMatthew G. Knepley } 5274b9ad928SBarry Smith PetscFunctionReturn(0); 5284b9ad928SBarry Smith } 5294b9ad928SBarry Smith 5306849ba73SBarry Smith static PetscErrorCode PCApplyTranspose_ASM(PC pc,Vec x,Vec y) 5314b9ad928SBarry Smith { 5324b9ad928SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 5336849ba73SBarry Smith PetscErrorCode ierr; 5341dd8081eSeaulisa PetscInt i,n_local_true = osm->n_local_true; 5354b9ad928SBarry Smith ScatterMode forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE; 5364b9ad928SBarry Smith 5374b9ad928SBarry Smith PetscFunctionBegin; 5384b9ad928SBarry Smith /* 5394b9ad928SBarry Smith Support for limiting the restriction or interpolation to only local 5404b9ad928SBarry Smith subdomain values (leaving the other values 0). 5414b9ad928SBarry Smith 5424b9ad928SBarry Smith Note: these are reversed from the PCApply_ASM() because we are applying the 5434b9ad928SBarry Smith transpose of the three terms 5444b9ad928SBarry Smith */ 545d8b3b5e3Seaulisa 5464b9ad928SBarry Smith if (!(osm->type & PC_ASM_INTERPOLATE)) { 5474b9ad928SBarry Smith forward = SCATTER_FORWARD_LOCAL; 5484b9ad928SBarry Smith /* have to zero the work RHS since scatter may leave some slots empty */ 5491dd8081eSeaulisa ierr = VecSet(osm->lx, 0.0);CHKERRQ(ierr); 5504b9ad928SBarry Smith } 5512fa5cd67SKarl Rupp if (!(osm->type & PC_ASM_RESTRICT)) reverse = SCATTER_REVERSE_LOCAL; 5524b9ad928SBarry Smith 553b0de9f23SBarry Smith /* zero the global and the local solutions */ 55410bd88b9SJed Brown ierr = VecZeroEntries(y);CHKERRQ(ierr); 555d8b3b5e3Seaulisa ierr = VecSet(osm->ly, 0.0);CHKERRQ(ierr); 556d8b3b5e3Seaulisa 557b0de9f23SBarry Smith /* Copy the global RHS to local RHS including the ghost nodes */ 5581dd8081eSeaulisa ierr = VecScatterBegin(osm->restriction, x, osm->lx, INSERT_VALUES, forward);CHKERRQ(ierr); 5591dd8081eSeaulisa ierr = VecScatterEnd(osm->restriction, x, osm->lx, INSERT_VALUES, forward);CHKERRQ(ierr); 560d8b3b5e3Seaulisa 561b0de9f23SBarry Smith /* Restrict local RHS to the overlapping 0-block RHS */ 5621dd8081eSeaulisa ierr = VecScatterBegin(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward);CHKERRQ(ierr); 5631dd8081eSeaulisa ierr = VecScatterEnd(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward);CHKERRQ(ierr); 564d8b3b5e3Seaulisa 5654b9ad928SBarry Smith /* do the local solves */ 566d8b3b5e3Seaulisa for (i = 0; i < n_local_true; ++i) { 567d8b3b5e3Seaulisa 568b0de9f23SBarry Smith /* solve the overlapping i-block */ 569d8b3b5e3Seaulisa ierr = KSPSolve(osm->ksp[i], osm->x[i], osm->y[i]);CHKERRQ(ierr); 570d8b3b5e3Seaulisa 571b0de9f23SBarry Smith if (osm->lprolongation) { /* interpolate the non-overalapping i-block solution to the local solution */ 5721dd8081eSeaulisa ierr = VecScatterBegin(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward);CHKERRQ(ierr); 5731dd8081eSeaulisa ierr = VecScatterEnd(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward);CHKERRQ(ierr); 574b9845e0eSMatthew Knepley } 575b0de9f23SBarry Smith else{ /* interpolate the overalapping i-block solution to the local solution */ 5761dd8081eSeaulisa ierr = VecScatterBegin(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse);CHKERRQ(ierr); 5771dd8081eSeaulisa ierr = VecScatterEnd(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse);CHKERRQ(ierr); 5784b9ad928SBarry Smith } 579d8b3b5e3Seaulisa 580d8b3b5e3Seaulisa if (i < n_local_true-1) { 581b0de9f23SBarry Smith /* Restrict local RHS to the overlapping (i+1)-block RHS */ 5821dd8081eSeaulisa ierr = VecScatterBegin(osm->lrestriction[i+1], osm->lx, osm->x[i+1], INSERT_VALUES, forward);CHKERRQ(ierr); 5831dd8081eSeaulisa ierr = VecScatterEnd(osm->lrestriction[i+1], osm->lx, osm->x[i+1], INSERT_VALUES, forward);CHKERRQ(ierr); 5844b9ad928SBarry Smith } 5854b9ad928SBarry Smith } 586b0de9f23SBarry Smith /* Add the local solution to the global solution including the ghost nodes */ 5871dd8081eSeaulisa ierr = VecScatterBegin(osm->restriction, osm->ly, y, ADD_VALUES, reverse);CHKERRQ(ierr); 5881dd8081eSeaulisa ierr = VecScatterEnd(osm->restriction, osm->ly, y, ADD_VALUES, reverse);CHKERRQ(ierr); 589d8b3b5e3Seaulisa 5904b9ad928SBarry Smith PetscFunctionReturn(0); 591d8b3b5e3Seaulisa 5924b9ad928SBarry Smith } 5934b9ad928SBarry Smith 594e91c6855SBarry Smith static PetscErrorCode PCReset_ASM(PC pc) 5954b9ad928SBarry Smith { 5964b9ad928SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 5976849ba73SBarry Smith PetscErrorCode ierr; 59813f74950SBarry Smith PetscInt i; 5994b9ad928SBarry Smith 6004b9ad928SBarry Smith PetscFunctionBegin; 60192bb6962SLisandro Dalcin if (osm->ksp) { 60292bb6962SLisandro Dalcin for (i=0; i<osm->n_local_true; i++) { 603e91c6855SBarry Smith ierr = KSPReset(osm->ksp[i]);CHKERRQ(ierr); 60492bb6962SLisandro Dalcin } 60592bb6962SLisandro Dalcin } 606e09e08ccSBarry Smith if (osm->pmat) { 60792bb6962SLisandro Dalcin if (osm->n_local_true > 0) { 60830a70a9aSHong Zhang ierr = MatDestroySubMatrices(osm->n_local_true,&osm->pmat);CHKERRQ(ierr); 60992bb6962SLisandro Dalcin } 61092bb6962SLisandro Dalcin } 6111dd8081eSeaulisa if (osm->lrestriction) { 6121dd8081eSeaulisa ierr = VecScatterDestroy(&osm->restriction);CHKERRQ(ierr); 6131dd8081eSeaulisa for (i=0; i<osm->n_local_true; i++) { 6141dd8081eSeaulisa ierr = VecScatterDestroy(&osm->lrestriction[i]);CHKERRQ(ierr); 6151dd8081eSeaulisa if (osm->lprolongation) {ierr = VecScatterDestroy(&osm->lprolongation[i]);CHKERRQ(ierr);} 616fcfd50ebSBarry Smith ierr = VecDestroy(&osm->x[i]);CHKERRQ(ierr); 617fcfd50ebSBarry Smith ierr = VecDestroy(&osm->y[i]);CHKERRQ(ierr); 6184b9ad928SBarry Smith } 6191dd8081eSeaulisa ierr = PetscFree(osm->lrestriction);CHKERRQ(ierr); 6201dd8081eSeaulisa if (osm->lprolongation) {ierr = PetscFree(osm->lprolongation);CHKERRQ(ierr);} 62105b42c5fSBarry Smith ierr = PetscFree(osm->x);CHKERRQ(ierr); 62278904715SLisandro Dalcin ierr = PetscFree(osm->y);CHKERRQ(ierr); 6231dd8081eSeaulisa 62492bb6962SLisandro Dalcin } 6252b691e39SMatthew Knepley ierr = PCASMDestroySubdomains(osm->n_local_true,osm->is,osm->is_local);CHKERRQ(ierr); 626fb745f2cSMatthew G. Knepley ierr = ISDestroy(&osm->lis);CHKERRQ(ierr); 627fb745f2cSMatthew G. Knepley ierr = VecDestroy(&osm->lx);CHKERRQ(ierr); 628fb745f2cSMatthew G. Knepley ierr = VecDestroy(&osm->ly);CHKERRQ(ierr); 629347574c9Seaulisa if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE) { 630347574c9Seaulisa ierr = MatDestroyMatrices(osm->n_local_true, &osm->lmats);CHKERRQ(ierr); 631fb745f2cSMatthew G. Knepley } 6322fa5cd67SKarl Rupp 63380ec0b7dSPatrick Sanan ierr = PetscFree(osm->sub_mat_type);CHKERRQ(ierr); 63480ec0b7dSPatrick Sanan 635e91c6855SBarry Smith osm->is = 0; 636e91c6855SBarry Smith osm->is_local = 0; 637e91c6855SBarry Smith PetscFunctionReturn(0); 638e91c6855SBarry Smith } 639e91c6855SBarry Smith 640e91c6855SBarry Smith static PetscErrorCode PCDestroy_ASM(PC pc) 641e91c6855SBarry Smith { 642e91c6855SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 643e91c6855SBarry Smith PetscErrorCode ierr; 644e91c6855SBarry Smith PetscInt i; 645e91c6855SBarry Smith 646e91c6855SBarry Smith PetscFunctionBegin; 647e91c6855SBarry Smith ierr = PCReset_ASM(pc);CHKERRQ(ierr); 648e91c6855SBarry Smith if (osm->ksp) { 649e91c6855SBarry Smith for (i=0; i<osm->n_local_true; i++) { 650fcfd50ebSBarry Smith ierr = KSPDestroy(&osm->ksp[i]);CHKERRQ(ierr); 651e91c6855SBarry Smith } 652e91c6855SBarry Smith ierr = PetscFree(osm->ksp);CHKERRQ(ierr); 653e91c6855SBarry Smith } 654e91c6855SBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 6554b9ad928SBarry Smith PetscFunctionReturn(0); 6564b9ad928SBarry Smith } 6574b9ad928SBarry Smith 6584416b707SBarry Smith static PetscErrorCode PCSetFromOptions_ASM(PetscOptionItems *PetscOptionsObject,PC pc) 6594b9ad928SBarry Smith { 6604b9ad928SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 6616849ba73SBarry Smith PetscErrorCode ierr; 6629dcbbd2bSBarry Smith PetscInt blocks,ovl; 663ace3abfcSBarry Smith PetscBool symset,flg; 66492bb6962SLisandro Dalcin PCASMType asmtype; 66512cd4985SMatthew G. Knepley PCCompositeType loctype; 66680ec0b7dSPatrick Sanan char sub_mat_type[256]; 6674b9ad928SBarry Smith 6684b9ad928SBarry Smith PetscFunctionBegin; 669bf108f30SBarry Smith /* set the type to symmetric if matrix is symmetric */ 67092bb6962SLisandro Dalcin if (!osm->type_set && pc->pmat) { 67192bb6962SLisandro Dalcin ierr = MatIsSymmetricKnown(pc->pmat,&symset,&flg);CHKERRQ(ierr); 6722fa5cd67SKarl Rupp if (symset && flg) osm->type = PC_ASM_BASIC; 673bf108f30SBarry Smith } 674e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"Additive Schwarz options");CHKERRQ(ierr); 675d709ea83SDmitry Karpeev ierr = PetscOptionsBool("-pc_asm_dm_subdomains","Use DMCreateDomainDecomposition() to define subdomains","PCASMSetDMSubdomains",osm->dm_subdomains,&osm->dm_subdomains,&flg);CHKERRQ(ierr); 67690d69ab7SBarry Smith ierr = PetscOptionsInt("-pc_asm_blocks","Number of subdomains","PCASMSetTotalSubdomains",osm->n,&blocks,&flg);CHKERRQ(ierr); 67765db9045SDmitry Karpeev if (flg) { 678f77a5242SKarl Rupp ierr = PCASMSetTotalSubdomains(pc,blocks,NULL,NULL);CHKERRQ(ierr); 679d709ea83SDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 68065db9045SDmitry Karpeev } 681*342c94f9SMatthew G. Knepley ierr = PetscOptionsInt("-pc_asm_local_blocks","Number of local subdomains","PCASMSetLocalSubdomains",osm->n_local_true,&blocks,&flg);CHKERRQ(ierr); 682*342c94f9SMatthew G. Knepley if (flg) { 683*342c94f9SMatthew G. Knepley ierr = PCASMSetLocalSubdomains(pc,blocks,NULL,NULL);CHKERRQ(ierr); 684*342c94f9SMatthew G. Knepley osm->dm_subdomains = PETSC_FALSE; 685*342c94f9SMatthew G. Knepley } 68690d69ab7SBarry Smith ierr = PetscOptionsInt("-pc_asm_overlap","Number of grid points overlap","PCASMSetOverlap",osm->overlap,&ovl,&flg);CHKERRQ(ierr); 68765db9045SDmitry Karpeev if (flg) { 68865db9045SDmitry Karpeev ierr = PCASMSetOverlap(pc,ovl);CHKERRQ(ierr); 689d709ea83SDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 69065db9045SDmitry Karpeev } 69190d69ab7SBarry Smith flg = PETSC_FALSE; 69290d69ab7SBarry Smith ierr = PetscOptionsEnum("-pc_asm_type","Type of restriction/extension","PCASMSetType",PCASMTypes,(PetscEnum)osm->type,(PetscEnum*)&asmtype,&flg);CHKERRQ(ierr); 69392bb6962SLisandro Dalcin if (flg) {ierr = PCASMSetType(pc,asmtype);CHKERRQ(ierr); } 69412cd4985SMatthew G. Knepley flg = PETSC_FALSE; 69512cd4985SMatthew G. Knepley ierr = PetscOptionsEnum("-pc_asm_local_type","Type of local solver composition","PCASMSetLocalType",PCCompositeTypes,(PetscEnum)osm->loctype,(PetscEnum*)&loctype,&flg);CHKERRQ(ierr); 69612cd4985SMatthew G. Knepley if (flg) {ierr = PCASMSetLocalType(pc,loctype);CHKERRQ(ierr); } 69780ec0b7dSPatrick Sanan ierr = PetscOptionsFList("-pc_asm_sub_mat_type","Subsolve Matrix Type","PCASMSetSubMatType",MatList,NULL,sub_mat_type,256,&flg);CHKERRQ(ierr); 69880ec0b7dSPatrick Sanan if(flg){ 699459726d8SSatish Balay ierr = PCASMSetSubMatType(pc,sub_mat_type);CHKERRQ(ierr); 70080ec0b7dSPatrick Sanan } 7014b9ad928SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 7024b9ad928SBarry Smith PetscFunctionReturn(0); 7034b9ad928SBarry Smith } 7044b9ad928SBarry Smith 7054b9ad928SBarry Smith /*------------------------------------------------------------------------------------*/ 7064b9ad928SBarry Smith 7071e6b0712SBarry Smith static PetscErrorCode PCASMSetLocalSubdomains_ASM(PC pc,PetscInt n,IS is[],IS is_local[]) 7084b9ad928SBarry Smith { 7094b9ad928SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 71092bb6962SLisandro Dalcin PetscErrorCode ierr; 71192bb6962SLisandro Dalcin PetscInt i; 7124b9ad928SBarry Smith 7134b9ad928SBarry Smith PetscFunctionBegin; 714e32f2f54SBarry Smith if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Each process must have 1 or more blocks, n = %D",n); 715ce94432eSBarry 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()."); 716e7e72b3dSBarry Smith 7174b9ad928SBarry Smith if (!pc->setupcalled) { 71892bb6962SLisandro Dalcin if (is) { 71992bb6962SLisandro Dalcin for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)is[i]);CHKERRQ(ierr);} 72092bb6962SLisandro Dalcin } 721832fc9a5SMatthew Knepley if (is_local) { 722832fc9a5SMatthew Knepley for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)is_local[i]);CHKERRQ(ierr);} 723832fc9a5SMatthew Knepley } 7242b691e39SMatthew Knepley ierr = PCASMDestroySubdomains(osm->n_local_true,osm->is,osm->is_local);CHKERRQ(ierr); 7252fa5cd67SKarl Rupp 7264b9ad928SBarry Smith osm->n_local_true = n; 72792bb6962SLisandro Dalcin osm->is = 0; 7282b691e39SMatthew Knepley osm->is_local = 0; 72992bb6962SLisandro Dalcin if (is) { 730785e854fSJed Brown ierr = PetscMalloc1(n,&osm->is);CHKERRQ(ierr); 7312fa5cd67SKarl Rupp for (i=0; i<n; i++) osm->is[i] = is[i]; 7323d03bb48SJed Brown /* Flag indicating that the user has set overlapping subdomains so PCASM should not increase their size. */ 7333d03bb48SJed Brown osm->overlap = -1; 73492bb6962SLisandro Dalcin } 7352b691e39SMatthew Knepley if (is_local) { 736785e854fSJed Brown ierr = PetscMalloc1(n,&osm->is_local);CHKERRQ(ierr); 7372fa5cd67SKarl Rupp for (i=0; i<n; i++) osm->is_local[i] = is_local[i]; 738a35b7b57SDmitry Karpeev if (!is) { 739785e854fSJed Brown ierr = PetscMalloc1(osm->n_local_true,&osm->is);CHKERRQ(ierr); 740a35b7b57SDmitry Karpeev for (i=0; i<osm->n_local_true; i++) { 741a35b7b57SDmitry Karpeev if (osm->overlap > 0) { /* With positive overlap, osm->is[i] will be modified */ 742a35b7b57SDmitry Karpeev ierr = ISDuplicate(osm->is_local[i],&osm->is[i]);CHKERRQ(ierr); 743a35b7b57SDmitry Karpeev ierr = ISCopy(osm->is_local[i],osm->is[i]);CHKERRQ(ierr); 744a35b7b57SDmitry Karpeev } else { 745a35b7b57SDmitry Karpeev ierr = PetscObjectReference((PetscObject)osm->is_local[i]);CHKERRQ(ierr); 746a35b7b57SDmitry Karpeev osm->is[i] = osm->is_local[i]; 747a35b7b57SDmitry Karpeev } 748a35b7b57SDmitry Karpeev } 749a35b7b57SDmitry Karpeev } 7502b691e39SMatthew Knepley } 7514b9ad928SBarry Smith } 7524b9ad928SBarry Smith PetscFunctionReturn(0); 7534b9ad928SBarry Smith } 7544b9ad928SBarry Smith 7551e6b0712SBarry Smith static PetscErrorCode PCASMSetTotalSubdomains_ASM(PC pc,PetscInt N,IS *is,IS *is_local) 7564b9ad928SBarry Smith { 7574b9ad928SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 7586849ba73SBarry Smith PetscErrorCode ierr; 75913f74950SBarry Smith PetscMPIInt rank,size; 76078904715SLisandro Dalcin PetscInt n; 7614b9ad928SBarry Smith 7624b9ad928SBarry Smith PetscFunctionBegin; 763ce94432eSBarry Smith if (N < 1) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Number of total blocks must be > 0, N = %D",N); 764ce94432eSBarry 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."); 7654b9ad928SBarry Smith 7664b9ad928SBarry Smith /* 767880770e9SJed Brown Split the subdomains equally among all processors 7684b9ad928SBarry Smith */ 769ce94432eSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);CHKERRQ(ierr); 770ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr); 7714b9ad928SBarry Smith n = N/size + ((N % size) > rank); 772e32f2f54SBarry 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); 773e32f2f54SBarry Smith if (pc->setupcalled && n != osm->n_local_true) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"PCASMSetTotalSubdomains() should be called before PCSetUp()."); 7744b9ad928SBarry Smith if (!pc->setupcalled) { 7752b691e39SMatthew Knepley ierr = PCASMDestroySubdomains(osm->n_local_true,osm->is,osm->is_local);CHKERRQ(ierr); 7762fa5cd67SKarl Rupp 7774b9ad928SBarry Smith osm->n_local_true = n; 7784b9ad928SBarry Smith osm->is = 0; 7792b691e39SMatthew Knepley osm->is_local = 0; 7804b9ad928SBarry Smith } 7814b9ad928SBarry Smith PetscFunctionReturn(0); 7824b9ad928SBarry Smith } 7834b9ad928SBarry Smith 7841e6b0712SBarry Smith static PetscErrorCode PCASMSetOverlap_ASM(PC pc,PetscInt ovl) 7854b9ad928SBarry Smith { 78692bb6962SLisandro Dalcin PC_ASM *osm = (PC_ASM*)pc->data; 7874b9ad928SBarry Smith 7884b9ad928SBarry Smith PetscFunctionBegin; 789ce94432eSBarry Smith if (ovl < 0) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap value requested"); 790ce94432eSBarry Smith if (pc->setupcalled && ovl != osm->overlap) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCASMSetOverlap() should be called before PCSetUp()."); 7912fa5cd67SKarl Rupp if (!pc->setupcalled) osm->overlap = ovl; 7924b9ad928SBarry Smith PetscFunctionReturn(0); 7934b9ad928SBarry Smith } 7944b9ad928SBarry Smith 7951e6b0712SBarry Smith static PetscErrorCode PCASMSetType_ASM(PC pc,PCASMType type) 7964b9ad928SBarry Smith { 79792bb6962SLisandro Dalcin PC_ASM *osm = (PC_ASM*)pc->data; 7984b9ad928SBarry Smith 7994b9ad928SBarry Smith PetscFunctionBegin; 8004b9ad928SBarry Smith osm->type = type; 801bf108f30SBarry Smith osm->type_set = PETSC_TRUE; 8024b9ad928SBarry Smith PetscFunctionReturn(0); 8034b9ad928SBarry Smith } 8044b9ad928SBarry Smith 805c60c7ad4SBarry Smith static PetscErrorCode PCASMGetType_ASM(PC pc,PCASMType *type) 806c60c7ad4SBarry Smith { 807c60c7ad4SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 808c60c7ad4SBarry Smith 809c60c7ad4SBarry Smith PetscFunctionBegin; 810c60c7ad4SBarry Smith *type = osm->type; 811c60c7ad4SBarry Smith PetscFunctionReturn(0); 812c60c7ad4SBarry Smith } 813c60c7ad4SBarry Smith 81412cd4985SMatthew G. Knepley static PetscErrorCode PCASMSetLocalType_ASM(PC pc, PCCompositeType type) 81512cd4985SMatthew G. Knepley { 81612cd4985SMatthew G. Knepley PC_ASM *osm = (PC_ASM *) pc->data; 81712cd4985SMatthew G. Knepley 81812cd4985SMatthew G. Knepley PetscFunctionBegin; 819d0ecd4dfSBarry 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"); 82012cd4985SMatthew G. Knepley osm->loctype = type; 82112cd4985SMatthew G. Knepley PetscFunctionReturn(0); 82212cd4985SMatthew G. Knepley } 82312cd4985SMatthew G. Knepley 82412cd4985SMatthew G. Knepley static PetscErrorCode PCASMGetLocalType_ASM(PC pc, PCCompositeType *type) 82512cd4985SMatthew G. Knepley { 82612cd4985SMatthew G. Knepley PC_ASM *osm = (PC_ASM *) pc->data; 82712cd4985SMatthew G. Knepley 82812cd4985SMatthew G. Knepley PetscFunctionBegin; 82912cd4985SMatthew G. Knepley *type = osm->loctype; 83012cd4985SMatthew G. Knepley PetscFunctionReturn(0); 83112cd4985SMatthew G. Knepley } 83212cd4985SMatthew G. Knepley 8331e6b0712SBarry Smith static PetscErrorCode PCASMSetSortIndices_ASM(PC pc,PetscBool doSort) 8346ed231c7SMatthew Knepley { 8356ed231c7SMatthew Knepley PC_ASM *osm = (PC_ASM*)pc->data; 8366ed231c7SMatthew Knepley 8376ed231c7SMatthew Knepley PetscFunctionBegin; 8386ed231c7SMatthew Knepley osm->sort_indices = doSort; 8396ed231c7SMatthew Knepley PetscFunctionReturn(0); 8406ed231c7SMatthew Knepley } 8416ed231c7SMatthew Knepley 8421e6b0712SBarry Smith static PetscErrorCode PCASMGetSubKSP_ASM(PC pc,PetscInt *n_local,PetscInt *first_local,KSP **ksp) 8434b9ad928SBarry Smith { 84492bb6962SLisandro Dalcin PC_ASM *osm = (PC_ASM*)pc->data; 845dfbe8321SBarry Smith PetscErrorCode ierr; 8464b9ad928SBarry Smith 8474b9ad928SBarry Smith PetscFunctionBegin; 848ce94432eSBarry 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"); 8494b9ad928SBarry Smith 8502fa5cd67SKarl Rupp if (n_local) *n_local = osm->n_local_true; 85192bb6962SLisandro Dalcin if (first_local) { 852ce94432eSBarry Smith ierr = MPI_Scan(&osm->n_local_true,first_local,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr); 85392bb6962SLisandro Dalcin *first_local -= osm->n_local_true; 85492bb6962SLisandro Dalcin } 85592bb6962SLisandro Dalcin if (ksp) { 85692bb6962SLisandro Dalcin /* Assume that local solves are now different; not necessarily 85792bb6962SLisandro Dalcin true though! This flag is used only for PCView_ASM() */ 85892bb6962SLisandro Dalcin *ksp = osm->ksp; 85992bb6962SLisandro Dalcin osm->same_local_solves = PETSC_FALSE; 86092bb6962SLisandro Dalcin } 8614b9ad928SBarry Smith PetscFunctionReturn(0); 8624b9ad928SBarry Smith } 8634b9ad928SBarry Smith 86480ec0b7dSPatrick Sanan static PetscErrorCode PCASMGetSubMatType_ASM(PC pc,MatType *sub_mat_type) 86580ec0b7dSPatrick Sanan { 86680ec0b7dSPatrick Sanan PC_ASM *osm = (PC_ASM*)pc->data; 86780ec0b7dSPatrick Sanan 86880ec0b7dSPatrick Sanan PetscFunctionBegin; 86980ec0b7dSPatrick Sanan PetscValidHeaderSpecific(pc,PC_CLASSID,1); 87080ec0b7dSPatrick Sanan PetscValidPointer(sub_mat_type,2); 87180ec0b7dSPatrick Sanan *sub_mat_type = osm->sub_mat_type; 87280ec0b7dSPatrick Sanan PetscFunctionReturn(0); 87380ec0b7dSPatrick Sanan } 87480ec0b7dSPatrick Sanan 87580ec0b7dSPatrick Sanan static PetscErrorCode PCASMSetSubMatType_ASM(PC pc,MatType sub_mat_type) 87680ec0b7dSPatrick Sanan { 87780ec0b7dSPatrick Sanan PetscErrorCode ierr; 87880ec0b7dSPatrick Sanan PC_ASM *osm = (PC_ASM*)pc->data; 87980ec0b7dSPatrick Sanan 88080ec0b7dSPatrick Sanan PetscFunctionBegin; 88180ec0b7dSPatrick Sanan PetscValidHeaderSpecific(pc,PC_CLASSID,1); 88280ec0b7dSPatrick Sanan ierr = PetscFree(osm->sub_mat_type);CHKERRQ(ierr); 88380ec0b7dSPatrick Sanan ierr = PetscStrallocpy(sub_mat_type,(char**)&osm->sub_mat_type);CHKERRQ(ierr); 88480ec0b7dSPatrick Sanan PetscFunctionReturn(0); 88580ec0b7dSPatrick Sanan } 88680ec0b7dSPatrick Sanan 8874b9ad928SBarry Smith /*@C 8881093a601SBarry Smith PCASMSetLocalSubdomains - Sets the local subdomains (for this processor only) for the additive Schwarz preconditioner. 8894b9ad928SBarry Smith 8904b9ad928SBarry Smith Collective on PC 8914b9ad928SBarry Smith 8924b9ad928SBarry Smith Input Parameters: 8934b9ad928SBarry Smith + pc - the preconditioner context 8944b9ad928SBarry Smith . n - the number of subdomains for this processor (default value = 1) 8958c03b21aSDmitry Karpeev . is - the index set that defines the subdomains for this processor 8960298fd71SBarry Smith (or NULL for PETSc to determine subdomains) 897f1ee410cSBarry 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 898f1ee410cSBarry Smith (or NULL to not provide these) 8994b9ad928SBarry Smith 900*342c94f9SMatthew G. Knepley Options Database Key: 901*342c94f9SMatthew G. Knepley To set the total number of subdomain blocks rather than specify the 902*342c94f9SMatthew G. Knepley index sets, use the option 903*342c94f9SMatthew G. Knepley . -pc_asm_local_blocks <blks> - Sets local blocks 904*342c94f9SMatthew G. Knepley 9054b9ad928SBarry Smith Notes: 9061093a601SBarry Smith The IS numbering is in the parallel, global numbering of the vector for both is and is_local 9074b9ad928SBarry Smith 9084b9ad928SBarry Smith By default the ASM preconditioner uses 1 block per processor. 9094b9ad928SBarry Smith 9104b9ad928SBarry Smith Use PCASMSetTotalSubdomains() to set the subdomains for all processors. 9114b9ad928SBarry Smith 912f1ee410cSBarry Smith If is_local is provided and PCASMType is PC_ASM_RESTRICT then the solution only over the is_local region is interpolated 913f1ee410cSBarry Smith back to form the global solution (this is the standard restricted additive Schwarz method) 914f1ee410cSBarry Smith 915f1ee410cSBarry 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 916f1ee410cSBarry Smith no code to handle that case. 917f1ee410cSBarry Smith 9184b9ad928SBarry Smith Level: advanced 9194b9ad928SBarry Smith 9204b9ad928SBarry Smith .keywords: PC, ASM, set, local, subdomains, additive Schwarz 9214b9ad928SBarry Smith 9224b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(), 923f1ee410cSBarry Smith PCASMCreateSubdomains2D(), PCASMGetLocalSubdomains(), PCASMType, PCASMSetType() 9244b9ad928SBarry Smith @*/ 9257087cfbeSBarry Smith PetscErrorCode PCASMSetLocalSubdomains(PC pc,PetscInt n,IS is[],IS is_local[]) 9264b9ad928SBarry Smith { 9277bb14e67SBarry Smith PetscErrorCode ierr; 9284b9ad928SBarry Smith 9294b9ad928SBarry Smith PetscFunctionBegin; 9300700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 9317bb14e67SBarry Smith ierr = PetscTryMethod(pc,"PCASMSetLocalSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,n,is,is_local));CHKERRQ(ierr); 9324b9ad928SBarry Smith PetscFunctionReturn(0); 9334b9ad928SBarry Smith } 9344b9ad928SBarry Smith 9354b9ad928SBarry Smith /*@C 936feb221f8SDmitry Karpeev PCASMSetTotalSubdomains - Sets the subdomains for all processors for the 9374b9ad928SBarry Smith additive Schwarz preconditioner. Either all or no processors in the 9384b9ad928SBarry Smith PC communicator must call this routine, with the same index sets. 9394b9ad928SBarry Smith 9404b9ad928SBarry Smith Collective on PC 9414b9ad928SBarry Smith 9424b9ad928SBarry Smith Input Parameters: 9434b9ad928SBarry Smith + pc - the preconditioner context 944feb221f8SDmitry Karpeev . N - the number of subdomains for all processors 945feb221f8SDmitry Karpeev . is - the index sets that define the subdomains for all processors 946dfaa0c5fSBarry Smith (or NULL to ask PETSc to determine the subdomains) 9472b691e39SMatthew Knepley - is_local - the index sets that define the local part of the subdomains for this processor 948f1ee410cSBarry Smith (or NULL to not provide this information) 9494b9ad928SBarry Smith 9504b9ad928SBarry Smith Options Database Key: 9514b9ad928SBarry Smith To set the total number of subdomain blocks rather than specify the 9524b9ad928SBarry Smith index sets, use the option 9534b9ad928SBarry Smith . -pc_asm_blocks <blks> - Sets total blocks 9544b9ad928SBarry Smith 9554b9ad928SBarry Smith Notes: 956f1ee410cSBarry Smith Currently you cannot use this to set the actual subdomains with the argument is or is_local. 9574b9ad928SBarry Smith 9584b9ad928SBarry Smith By default the ASM preconditioner uses 1 block per processor. 9594b9ad928SBarry Smith 9604b9ad928SBarry Smith These index sets cannot be destroyed until after completion of the 9614b9ad928SBarry Smith linear solves for which the ASM preconditioner is being used. 9624b9ad928SBarry Smith 9634b9ad928SBarry Smith Use PCASMSetLocalSubdomains() to set local subdomains. 9644b9ad928SBarry Smith 9651093a601SBarry Smith The IS numbering is in the parallel, global numbering of the vector for both is and is_local 9661093a601SBarry Smith 9674b9ad928SBarry Smith Level: advanced 9684b9ad928SBarry Smith 9694b9ad928SBarry Smith .keywords: PC, ASM, set, total, global, subdomains, additive Schwarz 9704b9ad928SBarry Smith 9714b9ad928SBarry Smith .seealso: PCASMSetLocalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(), 9724b9ad928SBarry Smith PCASMCreateSubdomains2D() 9734b9ad928SBarry Smith @*/ 9747087cfbeSBarry Smith PetscErrorCode PCASMSetTotalSubdomains(PC pc,PetscInt N,IS is[],IS is_local[]) 9754b9ad928SBarry Smith { 9767bb14e67SBarry Smith PetscErrorCode ierr; 9774b9ad928SBarry Smith 9784b9ad928SBarry Smith PetscFunctionBegin; 9790700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 9807bb14e67SBarry Smith ierr = PetscTryMethod(pc,"PCASMSetTotalSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,N,is,is_local));CHKERRQ(ierr); 9814b9ad928SBarry Smith PetscFunctionReturn(0); 9824b9ad928SBarry Smith } 9834b9ad928SBarry Smith 9844b9ad928SBarry Smith /*@ 9854b9ad928SBarry Smith PCASMSetOverlap - Sets the overlap between a pair of subdomains for the 9864b9ad928SBarry Smith additive Schwarz preconditioner. Either all or no processors in the 987f1ee410cSBarry Smith PC communicator must call this routine. 9884b9ad928SBarry Smith 989ad4df100SBarry Smith Logically Collective on PC 9904b9ad928SBarry Smith 9914b9ad928SBarry Smith Input Parameters: 9924b9ad928SBarry Smith + pc - the preconditioner context 9934b9ad928SBarry Smith - ovl - the amount of overlap between subdomains (ovl >= 0, default value = 1) 9944b9ad928SBarry Smith 9954b9ad928SBarry Smith Options Database Key: 9964b9ad928SBarry Smith . -pc_asm_overlap <ovl> - Sets overlap 9974b9ad928SBarry Smith 9984b9ad928SBarry Smith Notes: 9994b9ad928SBarry Smith By default the ASM preconditioner uses 1 block per processor. To use 10004b9ad928SBarry Smith multiple blocks per perocessor, see PCASMSetTotalSubdomains() and 10014b9ad928SBarry Smith PCASMSetLocalSubdomains() (and the option -pc_asm_blocks <blks>). 10024b9ad928SBarry Smith 10034b9ad928SBarry Smith The overlap defaults to 1, so if one desires that no additional 10044b9ad928SBarry Smith overlap be computed beyond what may have been set with a call to 10054b9ad928SBarry Smith PCASMSetTotalSubdomains() or PCASMSetLocalSubdomains(), then ovl 10064b9ad928SBarry Smith must be set to be 0. In particular, if one does not explicitly set 10074b9ad928SBarry Smith the subdomains an application code, then all overlap would be computed 10084b9ad928SBarry Smith internally by PETSc, and using an overlap of 0 would result in an ASM 10094b9ad928SBarry Smith variant that is equivalent to the block Jacobi preconditioner. 10104b9ad928SBarry Smith 1011f1ee410cSBarry Smith The default algorithm used by PETSc to increase overlap is fast, but not scalable, 1012f1ee410cSBarry Smith use the option -mat_increase_overlap_scalable when the problem and number of processes is large. 1013f1ee410cSBarry Smith 10144b9ad928SBarry Smith Note that one can define initial index sets with any overlap via 1015f1ee410cSBarry Smith PCASMSetLocalSubdomains(); the routine 10164b9ad928SBarry Smith PCASMSetOverlap() merely allows PETSc to extend that overlap further 10174b9ad928SBarry Smith if desired. 10184b9ad928SBarry Smith 10194b9ad928SBarry Smith Level: intermediate 10204b9ad928SBarry Smith 10214b9ad928SBarry Smith .keywords: PC, ASM, set, overlap 10224b9ad928SBarry Smith 10234b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetLocalSubdomains(), PCASMGetSubKSP(), 1024f1ee410cSBarry Smith PCASMCreateSubdomains2D(), PCASMGetLocalSubdomains(), MatIncreaseOverlap() 10254b9ad928SBarry Smith @*/ 10267087cfbeSBarry Smith PetscErrorCode PCASMSetOverlap(PC pc,PetscInt ovl) 10274b9ad928SBarry Smith { 10287bb14e67SBarry Smith PetscErrorCode ierr; 10294b9ad928SBarry Smith 10304b9ad928SBarry Smith PetscFunctionBegin; 10310700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1032c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,ovl,2); 10337bb14e67SBarry Smith ierr = PetscTryMethod(pc,"PCASMSetOverlap_C",(PC,PetscInt),(pc,ovl));CHKERRQ(ierr); 10344b9ad928SBarry Smith PetscFunctionReturn(0); 10354b9ad928SBarry Smith } 10364b9ad928SBarry Smith 10374b9ad928SBarry Smith /*@ 10384b9ad928SBarry Smith PCASMSetType - Sets the type of restriction and interpolation used 10394b9ad928SBarry Smith for local problems in the additive Schwarz method. 10404b9ad928SBarry Smith 1041ad4df100SBarry Smith Logically Collective on PC 10424b9ad928SBarry Smith 10434b9ad928SBarry Smith Input Parameters: 10444b9ad928SBarry Smith + pc - the preconditioner context 10454b9ad928SBarry Smith - type - variant of ASM, one of 10464b9ad928SBarry Smith .vb 10474b9ad928SBarry Smith PC_ASM_BASIC - full interpolation and restriction 10484b9ad928SBarry Smith PC_ASM_RESTRICT - full restriction, local processor interpolation 10494b9ad928SBarry Smith PC_ASM_INTERPOLATE - full interpolation, local processor restriction 10504b9ad928SBarry Smith PC_ASM_NONE - local processor restriction and interpolation 10514b9ad928SBarry Smith .ve 10524b9ad928SBarry Smith 10534b9ad928SBarry Smith Options Database Key: 10544b9ad928SBarry Smith . -pc_asm_type [basic,restrict,interpolate,none] - Sets ASM type 10554b9ad928SBarry Smith 1056f1ee410cSBarry Smith Notes: if the is_local arguments are passed to PCASMSetLocalSubdomains() then they are used when PC_ASM_RESTRICT has been selected 1057f1ee410cSBarry Smith to limit the local processor interpolation 1058f1ee410cSBarry Smith 10594b9ad928SBarry Smith Level: intermediate 10604b9ad928SBarry Smith 10614b9ad928SBarry Smith .keywords: PC, ASM, set, type 10624b9ad928SBarry Smith 10634b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetTotalSubdomains(), PCASMGetSubKSP(), 1064f1ee410cSBarry Smith PCASMCreateSubdomains2D(), PCASMType, PCASMSetLocalType(), PCASMGetLocalType() 10654b9ad928SBarry Smith @*/ 10667087cfbeSBarry Smith PetscErrorCode PCASMSetType(PC pc,PCASMType type) 10674b9ad928SBarry Smith { 10687bb14e67SBarry Smith PetscErrorCode ierr; 10694b9ad928SBarry Smith 10704b9ad928SBarry Smith PetscFunctionBegin; 10710700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1072c5eb9154SBarry Smith PetscValidLogicalCollectiveEnum(pc,type,2); 10737bb14e67SBarry Smith ierr = PetscTryMethod(pc,"PCASMSetType_C",(PC,PCASMType),(pc,type));CHKERRQ(ierr); 10744b9ad928SBarry Smith PetscFunctionReturn(0); 10754b9ad928SBarry Smith } 10764b9ad928SBarry Smith 1077c60c7ad4SBarry Smith /*@ 1078c60c7ad4SBarry Smith PCASMGetType - Gets the type of restriction and interpolation used 1079c60c7ad4SBarry Smith for local problems in the additive Schwarz method. 1080c60c7ad4SBarry Smith 1081c60c7ad4SBarry Smith Logically Collective on PC 1082c60c7ad4SBarry Smith 1083c60c7ad4SBarry Smith Input Parameter: 1084c60c7ad4SBarry Smith . pc - the preconditioner context 1085c60c7ad4SBarry Smith 1086c60c7ad4SBarry Smith Output Parameter: 1087c60c7ad4SBarry Smith . type - variant of ASM, one of 1088c60c7ad4SBarry Smith 1089c60c7ad4SBarry Smith .vb 1090c60c7ad4SBarry Smith PC_ASM_BASIC - full interpolation and restriction 1091c60c7ad4SBarry Smith PC_ASM_RESTRICT - full restriction, local processor interpolation 1092c60c7ad4SBarry Smith PC_ASM_INTERPOLATE - full interpolation, local processor restriction 1093c60c7ad4SBarry Smith PC_ASM_NONE - local processor restriction and interpolation 1094c60c7ad4SBarry Smith .ve 1095c60c7ad4SBarry Smith 1096c60c7ad4SBarry Smith Options Database Key: 1097c60c7ad4SBarry Smith . -pc_asm_type [basic,restrict,interpolate,none] - Sets ASM type 1098c60c7ad4SBarry Smith 1099c60c7ad4SBarry Smith Level: intermediate 1100c60c7ad4SBarry Smith 1101c60c7ad4SBarry Smith .keywords: PC, ASM, set, type 1102c60c7ad4SBarry Smith 1103c60c7ad4SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetTotalSubdomains(), PCASMGetSubKSP(), 1104f1ee410cSBarry Smith PCASMCreateSubdomains2D(), PCASMType, PCASMSetType(), PCASMSetLocalType(), PCASMGetLocalType() 1105c60c7ad4SBarry Smith @*/ 1106c60c7ad4SBarry Smith PetscErrorCode PCASMGetType(PC pc,PCASMType *type) 1107c60c7ad4SBarry Smith { 1108c60c7ad4SBarry Smith PetscErrorCode ierr; 1109c60c7ad4SBarry Smith 1110c60c7ad4SBarry Smith PetscFunctionBegin; 1111c60c7ad4SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1112c60c7ad4SBarry Smith ierr = PetscUseMethod(pc,"PCASMGetType_C",(PC,PCASMType*),(pc,type));CHKERRQ(ierr); 1113c60c7ad4SBarry Smith PetscFunctionReturn(0); 1114c60c7ad4SBarry Smith } 1115c60c7ad4SBarry Smith 111612cd4985SMatthew G. Knepley /*@ 111712cd4985SMatthew G. Knepley PCASMSetLocalType - Sets the type of composition used for local problems in the additive Schwarz method. 111812cd4985SMatthew G. Knepley 111912cd4985SMatthew G. Knepley Logically Collective on PC 112012cd4985SMatthew G. Knepley 112112cd4985SMatthew G. Knepley Input Parameters: 112212cd4985SMatthew G. Knepley + pc - the preconditioner context 112312cd4985SMatthew G. Knepley - type - type of composition, one of 112412cd4985SMatthew G. Knepley .vb 112512cd4985SMatthew G. Knepley PC_COMPOSITE_ADDITIVE - local additive combination 112612cd4985SMatthew G. Knepley PC_COMPOSITE_MULTIPLICATIVE - local multiplicative combination 112712cd4985SMatthew G. Knepley .ve 112812cd4985SMatthew G. Knepley 112912cd4985SMatthew G. Knepley Options Database Key: 113012cd4985SMatthew G. Knepley . -pc_asm_local_type [additive,multiplicative] - Sets local solver composition type 113112cd4985SMatthew G. Knepley 113212cd4985SMatthew G. Knepley Level: intermediate 113312cd4985SMatthew G. Knepley 1134f1ee410cSBarry Smith .seealso: PCASMSetType(), PCASMGetType(), PCASMGetLocalType(), PCASM, PCASMType, PCASMSetType(), PCASMGetType(), PCCompositeType 113512cd4985SMatthew G. Knepley @*/ 113612cd4985SMatthew G. Knepley PetscErrorCode PCASMSetLocalType(PC pc, PCCompositeType type) 113712cd4985SMatthew G. Knepley { 113812cd4985SMatthew G. Knepley PetscErrorCode ierr; 113912cd4985SMatthew G. Knepley 114012cd4985SMatthew G. Knepley PetscFunctionBegin; 114112cd4985SMatthew G. Knepley PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 114212cd4985SMatthew G. Knepley PetscValidLogicalCollectiveEnum(pc, type, 2); 114312cd4985SMatthew G. Knepley ierr = PetscTryMethod(pc, "PCASMSetLocalType_C", (PC, PCCompositeType), (pc, type));CHKERRQ(ierr); 114412cd4985SMatthew G. Knepley PetscFunctionReturn(0); 114512cd4985SMatthew G. Knepley } 114612cd4985SMatthew G. Knepley 114712cd4985SMatthew G. Knepley /*@ 114812cd4985SMatthew G. Knepley PCASMGetLocalType - Gets the type of composition used for local problems in the additive Schwarz method. 114912cd4985SMatthew G. Knepley 115012cd4985SMatthew G. Knepley Logically Collective on PC 115112cd4985SMatthew G. Knepley 115212cd4985SMatthew G. Knepley Input Parameter: 115312cd4985SMatthew G. Knepley . pc - the preconditioner context 115412cd4985SMatthew G. Knepley 115512cd4985SMatthew G. Knepley Output Parameter: 115612cd4985SMatthew G. Knepley . type - type of composition, one of 115712cd4985SMatthew G. Knepley .vb 115812cd4985SMatthew G. Knepley PC_COMPOSITE_ADDITIVE - local additive combination 115912cd4985SMatthew G. Knepley PC_COMPOSITE_MULTIPLICATIVE - local multiplicative combination 116012cd4985SMatthew G. Knepley .ve 116112cd4985SMatthew G. Knepley 116212cd4985SMatthew G. Knepley Options Database Key: 116312cd4985SMatthew G. Knepley . -pc_asm_local_type [additive,multiplicative] - Sets local solver composition type 116412cd4985SMatthew G. Knepley 116512cd4985SMatthew G. Knepley Level: intermediate 116612cd4985SMatthew G. Knepley 1167f1ee410cSBarry Smith .seealso: PCASMSetType(), PCASMGetType(), PCASMSetLocalType(), PCASMCreate(), PCASMType, PCASMSetType(), PCASMGetType(), PCCompositeType 116812cd4985SMatthew G. Knepley @*/ 116912cd4985SMatthew G. Knepley PetscErrorCode PCASMGetLocalType(PC pc, PCCompositeType *type) 117012cd4985SMatthew G. Knepley { 117112cd4985SMatthew G. Knepley PetscErrorCode ierr; 117212cd4985SMatthew G. Knepley 117312cd4985SMatthew G. Knepley PetscFunctionBegin; 117412cd4985SMatthew G. Knepley PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 117512cd4985SMatthew G. Knepley PetscValidPointer(type, 2); 117612cd4985SMatthew G. Knepley ierr = PetscUseMethod(pc, "PCASMGetLocalType_C", (PC, PCCompositeType *), (pc, type));CHKERRQ(ierr); 117712cd4985SMatthew G. Knepley PetscFunctionReturn(0); 117812cd4985SMatthew G. Knepley } 117912cd4985SMatthew G. Knepley 11806ed231c7SMatthew Knepley /*@ 11816ed231c7SMatthew Knepley PCASMSetSortIndices - Determines whether subdomain indices are sorted. 11826ed231c7SMatthew Knepley 1183ad4df100SBarry Smith Logically Collective on PC 11846ed231c7SMatthew Knepley 11856ed231c7SMatthew Knepley Input Parameters: 11866ed231c7SMatthew Knepley + pc - the preconditioner context 11876ed231c7SMatthew Knepley - doSort - sort the subdomain indices 11886ed231c7SMatthew Knepley 11896ed231c7SMatthew Knepley Level: intermediate 11906ed231c7SMatthew Knepley 11916ed231c7SMatthew Knepley .keywords: PC, ASM, set, type 11926ed231c7SMatthew Knepley 11936ed231c7SMatthew Knepley .seealso: PCASMSetLocalSubdomains(), PCASMSetTotalSubdomains(), PCASMGetSubKSP(), 11946ed231c7SMatthew Knepley PCASMCreateSubdomains2D() 11956ed231c7SMatthew Knepley @*/ 11967087cfbeSBarry Smith PetscErrorCode PCASMSetSortIndices(PC pc,PetscBool doSort) 11976ed231c7SMatthew Knepley { 11987bb14e67SBarry Smith PetscErrorCode ierr; 11996ed231c7SMatthew Knepley 12006ed231c7SMatthew Knepley PetscFunctionBegin; 12010700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1202acfcf0e5SJed Brown PetscValidLogicalCollectiveBool(pc,doSort,2); 12037bb14e67SBarry Smith ierr = PetscTryMethod(pc,"PCASMSetSortIndices_C",(PC,PetscBool),(pc,doSort));CHKERRQ(ierr); 12046ed231c7SMatthew Knepley PetscFunctionReturn(0); 12056ed231c7SMatthew Knepley } 12066ed231c7SMatthew Knepley 12074b9ad928SBarry Smith /*@C 12084b9ad928SBarry Smith PCASMGetSubKSP - Gets the local KSP contexts for all blocks on 12094b9ad928SBarry Smith this processor. 12104b9ad928SBarry Smith 12114b9ad928SBarry Smith Collective on PC iff first_local is requested 12124b9ad928SBarry Smith 12134b9ad928SBarry Smith Input Parameter: 12144b9ad928SBarry Smith . pc - the preconditioner context 12154b9ad928SBarry Smith 12164b9ad928SBarry Smith Output Parameters: 12170298fd71SBarry Smith + n_local - the number of blocks on this processor or NULL 12180298fd71SBarry Smith . first_local - the global number of the first block on this processor or NULL, 12190298fd71SBarry Smith all processors must request or all must pass NULL 12204b9ad928SBarry Smith - ksp - the array of KSP contexts 12214b9ad928SBarry Smith 12224b9ad928SBarry Smith Note: 1223d29017ddSJed Brown After PCASMGetSubKSP() the array of KSPes is not to be freed. 12244b9ad928SBarry Smith 12254b9ad928SBarry Smith You must call KSPSetUp() before calling PCASMGetSubKSP(). 12264b9ad928SBarry Smith 1227d29017ddSJed Brown Fortran note: 12282bf68e3eSBarry 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. 1229d29017ddSJed Brown 12304b9ad928SBarry Smith Level: advanced 12314b9ad928SBarry Smith 12324b9ad928SBarry Smith .keywords: PC, ASM, additive Schwarz, get, sub, KSP, context 12334b9ad928SBarry Smith 12344b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetTotalSubdomains(), PCASMSetOverlap(), 12354b9ad928SBarry Smith PCASMCreateSubdomains2D(), 12364b9ad928SBarry Smith @*/ 12377087cfbeSBarry Smith PetscErrorCode PCASMGetSubKSP(PC pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp[]) 12384b9ad928SBarry Smith { 12397bb14e67SBarry Smith PetscErrorCode ierr; 12404b9ad928SBarry Smith 12414b9ad928SBarry Smith PetscFunctionBegin; 12420700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 12437bb14e67SBarry Smith ierr = PetscUseMethod(pc,"PCASMGetSubKSP_C",(PC,PetscInt*,PetscInt*,KSP **),(pc,n_local,first_local,ksp));CHKERRQ(ierr); 12444b9ad928SBarry Smith PetscFunctionReturn(0); 12454b9ad928SBarry Smith } 12464b9ad928SBarry Smith 12474b9ad928SBarry Smith /* -------------------------------------------------------------------------------------*/ 12484b9ad928SBarry Smith /*MC 12493b09bd56SBarry Smith PCASM - Use the (restricted) additive Schwarz method, each block is (approximately) solved with 12504b9ad928SBarry Smith its own KSP object. 12514b9ad928SBarry Smith 12524b9ad928SBarry Smith Options Database Keys: 125349517cdeSBarry Smith + -pc_asm_blocks <blks> - Sets total blocks 12544b9ad928SBarry Smith . -pc_asm_overlap <ovl> - Sets overlap 1255f1ee410cSBarry Smith . -pc_asm_type [basic,restrict,interpolate,none] - Sets ASM type, default is restrict 1256f1ee410cSBarry Smith - -pc_asm_local_type [additive, multiplicative] - Sets ASM type, default is additive 12574b9ad928SBarry Smith 12583b09bd56SBarry Smith IMPORTANT: If you run with, for example, 3 blocks on 1 processor or 3 blocks on 3 processors you 12593b09bd56SBarry Smith will get a different convergence rate due to the default option of -pc_asm_type restrict. Use 12603b09bd56SBarry Smith -pc_asm_type basic to use the standard ASM. 12613b09bd56SBarry Smith 12624b9ad928SBarry Smith Notes: Each processor can have one or more blocks, but a block cannot be shared by more 1263f1ee410cSBarry Smith than one processor. Use PCGASM for subdomains shared by multiple processes. Defaults to one block per processor. 12644b9ad928SBarry Smith 12653b09bd56SBarry Smith To set options on the solvers for each block append -sub_ to all the KSP, and PC 1266d7ee0231SBarry Smith options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly 12674b9ad928SBarry Smith 1268a8c7a070SBarry Smith To set the options on the solvers separate for each block call PCASMGetSubKSP() 12694b9ad928SBarry Smith and set the options directly on the resulting KSP object (you can access its PC 12704b9ad928SBarry Smith with KSPGetPC()) 12714b9ad928SBarry Smith 12724b9ad928SBarry Smith Level: beginner 12734b9ad928SBarry Smith 12744b9ad928SBarry Smith Concepts: additive Schwarz method 12754b9ad928SBarry Smith 1276c582cd25SBarry Smith References: 127796a0c994SBarry Smith + 1. - M Dryja, OB Widlund, An additive variant of the Schwarz alternating method for the case of many subregions 127896a0c994SBarry Smith Courant Institute, New York University Technical report 127996a0c994SBarry Smith - 1. - Barry Smith, Petter Bjorstad, and William Gropp, Domain Decompositions: Parallel Multilevel Methods for Elliptic Partial Differential Equations, 128096a0c994SBarry Smith Cambridge University Press. 1281c582cd25SBarry Smith 12824b9ad928SBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 1283f1ee410cSBarry Smith PCBJACOBI, PCASMGetSubKSP(), PCASMSetLocalSubdomains(), PCASMType, PCASMGetType(), PCASMSetLocalType(), PCASMGetLocalType() 1284f1ee410cSBarry Smith PCASMSetTotalSubdomains(), PCSetModifySubmatrices(), PCASMSetOverlap(), PCASMSetType(), PCCompositeType 1285e09e08ccSBarry Smith 12864b9ad928SBarry Smith M*/ 12874b9ad928SBarry Smith 12888cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_ASM(PC pc) 12894b9ad928SBarry Smith { 1290dfbe8321SBarry Smith PetscErrorCode ierr; 12914b9ad928SBarry Smith PC_ASM *osm; 12924b9ad928SBarry Smith 12934b9ad928SBarry Smith PetscFunctionBegin; 1294b00a9115SJed Brown ierr = PetscNewLog(pc,&osm);CHKERRQ(ierr); 12952fa5cd67SKarl Rupp 12964b9ad928SBarry Smith osm->n = PETSC_DECIDE; 12974b9ad928SBarry Smith osm->n_local = 0; 12982b837212SDmitry Karpeev osm->n_local_true = PETSC_DECIDE; 12994b9ad928SBarry Smith osm->overlap = 1; 13004b9ad928SBarry Smith osm->ksp = 0; 13012b691e39SMatthew Knepley osm->restriction = 0; 13025b3c0d42Seaulisa osm->lprolongation = 0; 13031dd8081eSeaulisa osm->lrestriction = 0; 130492bb6962SLisandro Dalcin osm->x = 0; 130592bb6962SLisandro Dalcin osm->y = 0; 13064b9ad928SBarry Smith osm->is = 0; 13072b691e39SMatthew Knepley osm->is_local = 0; 13084b9ad928SBarry Smith osm->mat = 0; 13094b9ad928SBarry Smith osm->pmat = 0; 13104b9ad928SBarry Smith osm->type = PC_ASM_RESTRICT; 131112cd4985SMatthew G. Knepley osm->loctype = PC_COMPOSITE_ADDITIVE; 13124b9ad928SBarry Smith osm->same_local_solves = PETSC_TRUE; 13136ed231c7SMatthew Knepley osm->sort_indices = PETSC_TRUE; 1314d709ea83SDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 131580ec0b7dSPatrick Sanan osm->sub_mat_type = NULL; 13164b9ad928SBarry Smith 131792bb6962SLisandro Dalcin pc->data = (void*)osm; 13184b9ad928SBarry Smith pc->ops->apply = PCApply_ASM; 13194b9ad928SBarry Smith pc->ops->applytranspose = PCApplyTranspose_ASM; 13204b9ad928SBarry Smith pc->ops->setup = PCSetUp_ASM; 1321e91c6855SBarry Smith pc->ops->reset = PCReset_ASM; 13224b9ad928SBarry Smith pc->ops->destroy = PCDestroy_ASM; 13234b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_ASM; 13244b9ad928SBarry Smith pc->ops->setuponblocks = PCSetUpOnBlocks_ASM; 13254b9ad928SBarry Smith pc->ops->view = PCView_ASM; 13264b9ad928SBarry Smith pc->ops->applyrichardson = 0; 13274b9ad928SBarry Smith 1328bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetLocalSubdomains_C",PCASMSetLocalSubdomains_ASM);CHKERRQ(ierr); 1329bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetTotalSubdomains_C",PCASMSetTotalSubdomains_ASM);CHKERRQ(ierr); 1330bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetOverlap_C",PCASMSetOverlap_ASM);CHKERRQ(ierr); 1331bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetType_C",PCASMSetType_ASM);CHKERRQ(ierr); 1332c60c7ad4SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMGetType_C",PCASMGetType_ASM);CHKERRQ(ierr); 133312cd4985SMatthew G. Knepley ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetLocalType_C",PCASMSetLocalType_ASM);CHKERRQ(ierr); 133412cd4985SMatthew G. Knepley ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMGetLocalType_C",PCASMGetLocalType_ASM);CHKERRQ(ierr); 1335bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetSortIndices_C",PCASMSetSortIndices_ASM);CHKERRQ(ierr); 1336bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMGetSubKSP_C",PCASMGetSubKSP_ASM);CHKERRQ(ierr); 133780ec0b7dSPatrick Sanan ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMGetSubMatType_C",PCASMGetSubMatType_ASM);CHKERRQ(ierr); 133880ec0b7dSPatrick Sanan ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetSubMatType_C",PCASMSetSubMatType_ASM);CHKERRQ(ierr); 13394b9ad928SBarry Smith PetscFunctionReturn(0); 13404b9ad928SBarry Smith } 13414b9ad928SBarry Smith 134292bb6962SLisandro Dalcin /*@C 134392bb6962SLisandro Dalcin PCASMCreateSubdomains - Creates the index sets for the overlapping Schwarz 134492bb6962SLisandro Dalcin preconditioner for a any problem on a general grid. 134592bb6962SLisandro Dalcin 134692bb6962SLisandro Dalcin Collective 134792bb6962SLisandro Dalcin 134892bb6962SLisandro Dalcin Input Parameters: 134992bb6962SLisandro Dalcin + A - The global matrix operator 135092bb6962SLisandro Dalcin - n - the number of local blocks 135192bb6962SLisandro Dalcin 135292bb6962SLisandro Dalcin Output Parameters: 135392bb6962SLisandro Dalcin . outis - the array of index sets defining the subdomains 135492bb6962SLisandro Dalcin 135592bb6962SLisandro Dalcin Level: advanced 135692bb6962SLisandro Dalcin 13577d6bfa3bSBarry Smith Note: this generates nonoverlapping subdomains; the PCASM will generate the overlap 13587d6bfa3bSBarry Smith from these if you use PCASMSetLocalSubdomains() 13597d6bfa3bSBarry Smith 13607d6bfa3bSBarry Smith In the Fortran version you must provide the array outis[] already allocated of length n. 13617d6bfa3bSBarry Smith 136292bb6962SLisandro Dalcin .keywords: PC, ASM, additive Schwarz, create, subdomains, unstructured grid 136392bb6962SLisandro Dalcin 136492bb6962SLisandro Dalcin .seealso: PCASMSetLocalSubdomains(), PCASMDestroySubdomains() 136592bb6962SLisandro Dalcin @*/ 13667087cfbeSBarry Smith PetscErrorCode PCASMCreateSubdomains(Mat A, PetscInt n, IS* outis[]) 136792bb6962SLisandro Dalcin { 136892bb6962SLisandro Dalcin MatPartitioning mpart; 136992bb6962SLisandro Dalcin const char *prefix; 1370c56a70eeSBarry Smith void (*f)(void); 137192bb6962SLisandro Dalcin PetscInt i,j,rstart,rend,bs; 137211bd1e4dSLisandro Dalcin PetscBool isbaij = PETSC_FALSE,foundpart = PETSC_FALSE; 13730298fd71SBarry Smith Mat Ad = NULL, adj; 137492bb6962SLisandro Dalcin IS ispart,isnumb,*is; 137592bb6962SLisandro Dalcin PetscErrorCode ierr; 137692bb6962SLisandro Dalcin 137792bb6962SLisandro Dalcin PetscFunctionBegin; 13780700a824SBarry Smith PetscValidHeaderSpecific(A,MAT_CLASSID,1); 137992bb6962SLisandro Dalcin PetscValidPointer(outis,3); 1380c1235816SBarry Smith if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of local blocks must be > 0, n = %D",n); 138192bb6962SLisandro Dalcin 138292bb6962SLisandro Dalcin /* Get prefix, row distribution, and block size */ 138392bb6962SLisandro Dalcin ierr = MatGetOptionsPrefix(A,&prefix);CHKERRQ(ierr); 138492bb6962SLisandro Dalcin ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 138592bb6962SLisandro Dalcin ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr); 138665e19b50SBarry 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); 138765e19b50SBarry Smith 138892bb6962SLisandro Dalcin /* Get diagonal block from matrix if possible */ 1389c56a70eeSBarry Smith ierr = MatShellGetOperation(A,MATOP_GET_DIAGONAL_BLOCK,&f);CHKERRQ(ierr); 139092bb6962SLisandro Dalcin if (f) { 139111bd1e4dSLisandro Dalcin ierr = MatGetDiagonalBlock(A,&Ad);CHKERRQ(ierr); 139292bb6962SLisandro Dalcin } 139392bb6962SLisandro Dalcin if (Ad) { 1394251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)Ad,MATSEQBAIJ,&isbaij);CHKERRQ(ierr); 1395251f4c67SDmitry Karpeev if (!isbaij) {ierr = PetscObjectTypeCompare((PetscObject)Ad,MATSEQSBAIJ,&isbaij);CHKERRQ(ierr);} 139692bb6962SLisandro Dalcin } 139792bb6962SLisandro Dalcin if (Ad && n > 1) { 1398ace3abfcSBarry Smith PetscBool match,done; 139992bb6962SLisandro Dalcin /* Try to setup a good matrix partitioning if available */ 140092bb6962SLisandro Dalcin ierr = MatPartitioningCreate(PETSC_COMM_SELF,&mpart);CHKERRQ(ierr); 140192bb6962SLisandro Dalcin ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr); 140292bb6962SLisandro Dalcin ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr); 1403251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGCURRENT,&match);CHKERRQ(ierr); 140492bb6962SLisandro Dalcin if (!match) { 1405251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGSQUARE,&match);CHKERRQ(ierr); 140692bb6962SLisandro Dalcin } 140792bb6962SLisandro Dalcin if (!match) { /* assume a "good" partitioner is available */ 14081a83f524SJed Brown PetscInt na; 14091a83f524SJed Brown const PetscInt *ia,*ja; 141092bb6962SLisandro Dalcin ierr = MatGetRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr); 141192bb6962SLisandro Dalcin if (done) { 141292bb6962SLisandro Dalcin /* Build adjacency matrix by hand. Unfortunately a call to 141392bb6962SLisandro Dalcin MatConvert(Ad,MATMPIADJ,MAT_INITIAL_MATRIX,&adj) will 141492bb6962SLisandro Dalcin remove the block-aij structure and we cannot expect 141592bb6962SLisandro Dalcin MatPartitioning to split vertices as we need */ 14161a83f524SJed Brown PetscInt i,j,len,nnz,cnt,*iia=0,*jja=0; 14171a83f524SJed Brown const PetscInt *row; 141892bb6962SLisandro Dalcin nnz = 0; 141992bb6962SLisandro Dalcin for (i=0; i<na; i++) { /* count number of nonzeros */ 142092bb6962SLisandro Dalcin len = ia[i+1] - ia[i]; 142192bb6962SLisandro Dalcin row = ja + ia[i]; 142292bb6962SLisandro Dalcin for (j=0; j<len; j++) { 142392bb6962SLisandro Dalcin if (row[j] == i) { /* don't count diagonal */ 142492bb6962SLisandro Dalcin len--; break; 142592bb6962SLisandro Dalcin } 142692bb6962SLisandro Dalcin } 142792bb6962SLisandro Dalcin nnz += len; 142892bb6962SLisandro Dalcin } 1429854ce69bSBarry Smith ierr = PetscMalloc1(na+1,&iia);CHKERRQ(ierr); 1430854ce69bSBarry Smith ierr = PetscMalloc1(nnz,&jja);CHKERRQ(ierr); 143192bb6962SLisandro Dalcin nnz = 0; 143292bb6962SLisandro Dalcin iia[0] = 0; 143392bb6962SLisandro Dalcin for (i=0; i<na; i++) { /* fill adjacency */ 143492bb6962SLisandro Dalcin cnt = 0; 143592bb6962SLisandro Dalcin len = ia[i+1] - ia[i]; 143692bb6962SLisandro Dalcin row = ja + ia[i]; 143792bb6962SLisandro Dalcin for (j=0; j<len; j++) { 143892bb6962SLisandro Dalcin if (row[j] != i) { /* if not diagonal */ 143992bb6962SLisandro Dalcin jja[nnz+cnt++] = row[j]; 144092bb6962SLisandro Dalcin } 144192bb6962SLisandro Dalcin } 144292bb6962SLisandro Dalcin nnz += cnt; 144392bb6962SLisandro Dalcin iia[i+1] = nnz; 144492bb6962SLisandro Dalcin } 144592bb6962SLisandro Dalcin /* Partitioning of the adjacency matrix */ 14460298fd71SBarry Smith ierr = MatCreateMPIAdj(PETSC_COMM_SELF,na,na,iia,jja,NULL,&adj);CHKERRQ(ierr); 144792bb6962SLisandro Dalcin ierr = MatPartitioningSetAdjacency(mpart,adj);CHKERRQ(ierr); 144892bb6962SLisandro Dalcin ierr = MatPartitioningSetNParts(mpart,n);CHKERRQ(ierr); 144992bb6962SLisandro Dalcin ierr = MatPartitioningApply(mpart,&ispart);CHKERRQ(ierr); 145092bb6962SLisandro Dalcin ierr = ISPartitioningToNumbering(ispart,&isnumb);CHKERRQ(ierr); 1451fcfd50ebSBarry Smith ierr = MatDestroy(&adj);CHKERRQ(ierr); 145292bb6962SLisandro Dalcin foundpart = PETSC_TRUE; 145392bb6962SLisandro Dalcin } 145492bb6962SLisandro Dalcin ierr = MatRestoreRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr); 145592bb6962SLisandro Dalcin } 1456fcfd50ebSBarry Smith ierr = MatPartitioningDestroy(&mpart);CHKERRQ(ierr); 145792bb6962SLisandro Dalcin } 145892bb6962SLisandro Dalcin 1459785e854fSJed Brown ierr = PetscMalloc1(n,&is);CHKERRQ(ierr); 146092bb6962SLisandro Dalcin *outis = is; 146192bb6962SLisandro Dalcin 146292bb6962SLisandro Dalcin if (!foundpart) { 146392bb6962SLisandro Dalcin 146492bb6962SLisandro Dalcin /* Partitioning by contiguous chunks of rows */ 146592bb6962SLisandro Dalcin 146692bb6962SLisandro Dalcin PetscInt mbs = (rend-rstart)/bs; 146792bb6962SLisandro Dalcin PetscInt start = rstart; 146892bb6962SLisandro Dalcin for (i=0; i<n; i++) { 146992bb6962SLisandro Dalcin PetscInt count = (mbs/n + ((mbs % n) > i)) * bs; 147092bb6962SLisandro Dalcin ierr = ISCreateStride(PETSC_COMM_SELF,count,start,1,&is[i]);CHKERRQ(ierr); 147192bb6962SLisandro Dalcin start += count; 147292bb6962SLisandro Dalcin } 147392bb6962SLisandro Dalcin 147492bb6962SLisandro Dalcin } else { 147592bb6962SLisandro Dalcin 147692bb6962SLisandro Dalcin /* Partitioning by adjacency of diagonal block */ 147792bb6962SLisandro Dalcin 147892bb6962SLisandro Dalcin const PetscInt *numbering; 147992bb6962SLisandro Dalcin PetscInt *count,nidx,*indices,*newidx,start=0; 148092bb6962SLisandro Dalcin /* Get node count in each partition */ 1481785e854fSJed Brown ierr = PetscMalloc1(n,&count);CHKERRQ(ierr); 148292bb6962SLisandro Dalcin ierr = ISPartitioningCount(ispart,n,count);CHKERRQ(ierr); 148392bb6962SLisandro Dalcin if (isbaij && bs > 1) { /* adjust for the block-aij case */ 148492bb6962SLisandro Dalcin for (i=0; i<n; i++) count[i] *= bs; 148592bb6962SLisandro Dalcin } 148692bb6962SLisandro Dalcin /* Build indices from node numbering */ 148792bb6962SLisandro Dalcin ierr = ISGetLocalSize(isnumb,&nidx);CHKERRQ(ierr); 1488785e854fSJed Brown ierr = PetscMalloc1(nidx,&indices);CHKERRQ(ierr); 148992bb6962SLisandro Dalcin for (i=0; i<nidx; i++) indices[i] = i; /* needs to be initialized */ 149092bb6962SLisandro Dalcin ierr = ISGetIndices(isnumb,&numbering);CHKERRQ(ierr); 149192bb6962SLisandro Dalcin ierr = PetscSortIntWithPermutation(nidx,numbering,indices);CHKERRQ(ierr); 149292bb6962SLisandro Dalcin ierr = ISRestoreIndices(isnumb,&numbering);CHKERRQ(ierr); 149392bb6962SLisandro Dalcin if (isbaij && bs > 1) { /* adjust for the block-aij case */ 1494785e854fSJed Brown ierr = PetscMalloc1(nidx*bs,&newidx);CHKERRQ(ierr); 14952fa5cd67SKarl Rupp for (i=0; i<nidx; i++) { 14962fa5cd67SKarl Rupp for (j=0; j<bs; j++) newidx[i*bs+j] = indices[i]*bs + j; 14972fa5cd67SKarl Rupp } 149892bb6962SLisandro Dalcin ierr = PetscFree(indices);CHKERRQ(ierr); 149992bb6962SLisandro Dalcin nidx *= bs; 150092bb6962SLisandro Dalcin indices = newidx; 150192bb6962SLisandro Dalcin } 150292bb6962SLisandro Dalcin /* Shift to get global indices */ 150392bb6962SLisandro Dalcin for (i=0; i<nidx; i++) indices[i] += rstart; 150492bb6962SLisandro Dalcin 150592bb6962SLisandro Dalcin /* Build the index sets for each block */ 150692bb6962SLisandro Dalcin for (i=0; i<n; i++) { 150770b3c8c7SBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,count[i],&indices[start],PETSC_COPY_VALUES,&is[i]);CHKERRQ(ierr); 150892bb6962SLisandro Dalcin ierr = ISSort(is[i]);CHKERRQ(ierr); 150992bb6962SLisandro Dalcin start += count[i]; 151092bb6962SLisandro Dalcin } 151192bb6962SLisandro Dalcin 15123bf036e2SBarry Smith ierr = PetscFree(count);CHKERRQ(ierr); 15133bf036e2SBarry Smith ierr = PetscFree(indices);CHKERRQ(ierr); 1514fcfd50ebSBarry Smith ierr = ISDestroy(&isnumb);CHKERRQ(ierr); 1515fcfd50ebSBarry Smith ierr = ISDestroy(&ispart);CHKERRQ(ierr); 151692bb6962SLisandro Dalcin 151792bb6962SLisandro Dalcin } 151892bb6962SLisandro Dalcin PetscFunctionReturn(0); 151992bb6962SLisandro Dalcin } 152092bb6962SLisandro Dalcin 152192bb6962SLisandro Dalcin /*@C 152292bb6962SLisandro Dalcin PCASMDestroySubdomains - Destroys the index sets created with 152392bb6962SLisandro Dalcin PCASMCreateSubdomains(). Should be called after setting subdomains 152492bb6962SLisandro Dalcin with PCASMSetLocalSubdomains(). 152592bb6962SLisandro Dalcin 152692bb6962SLisandro Dalcin Collective 152792bb6962SLisandro Dalcin 152892bb6962SLisandro Dalcin Input Parameters: 152992bb6962SLisandro Dalcin + n - the number of index sets 15302b691e39SMatthew Knepley . is - the array of index sets 15310298fd71SBarry Smith - is_local - the array of local index sets, can be NULL 153292bb6962SLisandro Dalcin 153392bb6962SLisandro Dalcin Level: advanced 153492bb6962SLisandro Dalcin 153592bb6962SLisandro Dalcin .keywords: PC, ASM, additive Schwarz, create, subdomains, unstructured grid 153692bb6962SLisandro Dalcin 153792bb6962SLisandro Dalcin .seealso: PCASMCreateSubdomains(), PCASMSetLocalSubdomains() 153892bb6962SLisandro Dalcin @*/ 15397087cfbeSBarry Smith PetscErrorCode PCASMDestroySubdomains(PetscInt n, IS is[], IS is_local[]) 154092bb6962SLisandro Dalcin { 154192bb6962SLisandro Dalcin PetscInt i; 154292bb6962SLisandro Dalcin PetscErrorCode ierr; 15435fd66863SKarl Rupp 154492bb6962SLisandro Dalcin PetscFunctionBegin; 1545a35b7b57SDmitry Karpeev if (n <= 0) PetscFunctionReturn(0); 1546a35b7b57SDmitry Karpeev if (is) { 154792bb6962SLisandro Dalcin PetscValidPointer(is,2); 1548fcfd50ebSBarry Smith for (i=0; i<n; i++) { ierr = ISDestroy(&is[i]);CHKERRQ(ierr); } 154992bb6962SLisandro Dalcin ierr = PetscFree(is);CHKERRQ(ierr); 1550a35b7b57SDmitry Karpeev } 15512b691e39SMatthew Knepley if (is_local) { 15522b691e39SMatthew Knepley PetscValidPointer(is_local,3); 1553fcfd50ebSBarry Smith for (i=0; i<n; i++) { ierr = ISDestroy(&is_local[i]);CHKERRQ(ierr); } 15542b691e39SMatthew Knepley ierr = PetscFree(is_local);CHKERRQ(ierr); 15552b691e39SMatthew Knepley } 155692bb6962SLisandro Dalcin PetscFunctionReturn(0); 155792bb6962SLisandro Dalcin } 155892bb6962SLisandro Dalcin 15594b9ad928SBarry Smith /*@ 15604b9ad928SBarry Smith PCASMCreateSubdomains2D - Creates the index sets for the overlapping Schwarz 15614b9ad928SBarry Smith preconditioner for a two-dimensional problem on a regular grid. 15624b9ad928SBarry Smith 15634b9ad928SBarry Smith Not Collective 15644b9ad928SBarry Smith 15654b9ad928SBarry Smith Input Parameters: 15664b9ad928SBarry Smith + m, n - the number of mesh points in the x and y directions 15674b9ad928SBarry Smith . M, N - the number of subdomains in the x and y directions 15684b9ad928SBarry Smith . dof - degrees of freedom per node 15694b9ad928SBarry Smith - overlap - overlap in mesh lines 15704b9ad928SBarry Smith 15714b9ad928SBarry Smith Output Parameters: 15724b9ad928SBarry Smith + Nsub - the number of subdomains created 15733d03bb48SJed Brown . is - array of index sets defining overlapping (if overlap > 0) subdomains 15743d03bb48SJed Brown - is_local - array of index sets defining non-overlapping subdomains 15754b9ad928SBarry Smith 15764b9ad928SBarry Smith Note: 15774b9ad928SBarry Smith Presently PCAMSCreateSubdomains2d() is valid only for sequential 15784b9ad928SBarry Smith preconditioners. More general related routines are 15794b9ad928SBarry Smith PCASMSetTotalSubdomains() and PCASMSetLocalSubdomains(). 15804b9ad928SBarry Smith 15814b9ad928SBarry Smith Level: advanced 15824b9ad928SBarry Smith 15834b9ad928SBarry Smith .keywords: PC, ASM, additive Schwarz, create, subdomains, 2D, regular grid 15844b9ad928SBarry Smith 15854b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetLocalSubdomains(), PCASMGetSubKSP(), 15864b9ad928SBarry Smith PCASMSetOverlap() 15874b9ad928SBarry Smith @*/ 15887087cfbeSBarry Smith PetscErrorCode PCASMCreateSubdomains2D(PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt dof,PetscInt overlap,PetscInt *Nsub,IS **is,IS **is_local) 15894b9ad928SBarry Smith { 15903d03bb48SJed Brown PetscInt i,j,height,width,ystart,xstart,yleft,yright,xleft,xright,loc_outer; 15916849ba73SBarry Smith PetscErrorCode ierr; 159213f74950SBarry Smith PetscInt nidx,*idx,loc,ii,jj,count; 15934b9ad928SBarry Smith 15944b9ad928SBarry Smith PetscFunctionBegin; 1595e32f2f54SBarry Smith if (dof != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP," "); 15964b9ad928SBarry Smith 15974b9ad928SBarry Smith *Nsub = N*M; 1598854ce69bSBarry Smith ierr = PetscMalloc1(*Nsub,is);CHKERRQ(ierr); 1599854ce69bSBarry Smith ierr = PetscMalloc1(*Nsub,is_local);CHKERRQ(ierr); 16004b9ad928SBarry Smith ystart = 0; 16013d03bb48SJed Brown loc_outer = 0; 16024b9ad928SBarry Smith for (i=0; i<N; i++) { 16034b9ad928SBarry Smith height = n/N + ((n % N) > i); /* height of subdomain */ 1604e32f2f54SBarry Smith if (height < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many N subdomains for mesh dimension n"); 16054b9ad928SBarry Smith yleft = ystart - overlap; if (yleft < 0) yleft = 0; 16064b9ad928SBarry Smith yright = ystart + height + overlap; if (yright > n) yright = n; 16074b9ad928SBarry Smith xstart = 0; 16084b9ad928SBarry Smith for (j=0; j<M; j++) { 16094b9ad928SBarry Smith width = m/M + ((m % M) > j); /* width of subdomain */ 1610e32f2f54SBarry Smith if (width < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many M subdomains for mesh dimension m"); 16114b9ad928SBarry Smith xleft = xstart - overlap; if (xleft < 0) xleft = 0; 16124b9ad928SBarry Smith xright = xstart + width + overlap; if (xright > m) xright = m; 16134b9ad928SBarry Smith nidx = (xright - xleft)*(yright - yleft); 1614785e854fSJed Brown ierr = PetscMalloc1(nidx,&idx);CHKERRQ(ierr); 16154b9ad928SBarry Smith loc = 0; 16164b9ad928SBarry Smith for (ii=yleft; ii<yright; ii++) { 16174b9ad928SBarry Smith count = m*ii + xleft; 16182fa5cd67SKarl Rupp for (jj=xleft; jj<xright; jj++) idx[loc++] = count++; 16194b9ad928SBarry Smith } 162070b3c8c7SBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,nidx,idx,PETSC_COPY_VALUES,(*is)+loc_outer);CHKERRQ(ierr); 16213d03bb48SJed Brown if (overlap == 0) { 16223d03bb48SJed Brown ierr = PetscObjectReference((PetscObject)(*is)[loc_outer]);CHKERRQ(ierr); 16232fa5cd67SKarl Rupp 16243d03bb48SJed Brown (*is_local)[loc_outer] = (*is)[loc_outer]; 16253d03bb48SJed Brown } else { 16263d03bb48SJed Brown for (loc=0,ii=ystart; ii<ystart+height; ii++) { 16273d03bb48SJed Brown for (jj=xstart; jj<xstart+width; jj++) { 16283d03bb48SJed Brown idx[loc++] = m*ii + jj; 16293d03bb48SJed Brown } 16303d03bb48SJed Brown } 163170b3c8c7SBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,loc,idx,PETSC_COPY_VALUES,*is_local+loc_outer);CHKERRQ(ierr); 16323d03bb48SJed Brown } 16334b9ad928SBarry Smith ierr = PetscFree(idx);CHKERRQ(ierr); 16344b9ad928SBarry Smith xstart += width; 16353d03bb48SJed Brown loc_outer++; 16364b9ad928SBarry Smith } 16374b9ad928SBarry Smith ystart += height; 16384b9ad928SBarry Smith } 16394b9ad928SBarry Smith for (i=0; i<*Nsub; i++) { ierr = ISSort((*is)[i]);CHKERRQ(ierr); } 16404b9ad928SBarry Smith PetscFunctionReturn(0); 16414b9ad928SBarry Smith } 16424b9ad928SBarry Smith 16434b9ad928SBarry Smith /*@C 16444b9ad928SBarry Smith PCASMGetLocalSubdomains - Gets the local subdomains (for this processor 16454b9ad928SBarry Smith only) for the additive Schwarz preconditioner. 16464b9ad928SBarry Smith 1647ad4df100SBarry Smith Not Collective 16484b9ad928SBarry Smith 16494b9ad928SBarry Smith Input Parameter: 16504b9ad928SBarry Smith . pc - the preconditioner context 16514b9ad928SBarry Smith 16524b9ad928SBarry Smith Output Parameters: 16534b9ad928SBarry Smith + n - the number of subdomains for this processor (default value = 1) 16542b691e39SMatthew Knepley . is - the index sets that define the subdomains for this processor 16550298fd71SBarry Smith - is_local - the index sets that define the local part of the subdomains for this processor (can be NULL) 16564b9ad928SBarry Smith 16574b9ad928SBarry Smith 16584b9ad928SBarry Smith Notes: 16594b9ad928SBarry Smith The IS numbering is in the parallel, global numbering of the vector. 16604b9ad928SBarry Smith 16614b9ad928SBarry Smith Level: advanced 16624b9ad928SBarry Smith 16634b9ad928SBarry Smith .keywords: PC, ASM, set, local, subdomains, additive Schwarz 16644b9ad928SBarry Smith 16654b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(), 16664b9ad928SBarry Smith PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubmatrices() 16674b9ad928SBarry Smith @*/ 16687087cfbeSBarry Smith PetscErrorCode PCASMGetLocalSubdomains(PC pc,PetscInt *n,IS *is[],IS *is_local[]) 16694b9ad928SBarry Smith { 16702a808120SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 167192bb6962SLisandro Dalcin PetscErrorCode ierr; 1672ace3abfcSBarry Smith PetscBool match; 16734b9ad928SBarry Smith 16744b9ad928SBarry Smith PetscFunctionBegin; 16750700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 16764482741eSBarry Smith PetscValidIntPointer(n,2); 167792bb6962SLisandro Dalcin if (is) PetscValidPointer(is,3); 1678251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)pc,PCASM,&match);CHKERRQ(ierr); 16792a808120SBarry Smith if (!match) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"PC is not a PCASM"); 16804b9ad928SBarry Smith if (n) *n = osm->n_local_true; 16814b9ad928SBarry Smith if (is) *is = osm->is; 16822b691e39SMatthew Knepley if (is_local) *is_local = osm->is_local; 16834b9ad928SBarry Smith PetscFunctionReturn(0); 16844b9ad928SBarry Smith } 16854b9ad928SBarry Smith 16864b9ad928SBarry Smith /*@C 16874b9ad928SBarry Smith PCASMGetLocalSubmatrices - Gets the local submatrices (for this processor 16884b9ad928SBarry Smith only) for the additive Schwarz preconditioner. 16894b9ad928SBarry Smith 1690ad4df100SBarry Smith Not Collective 16914b9ad928SBarry Smith 16924b9ad928SBarry Smith Input Parameter: 16934b9ad928SBarry Smith . pc - the preconditioner context 16944b9ad928SBarry Smith 16954b9ad928SBarry Smith Output Parameters: 16964b9ad928SBarry Smith + n - the number of matrices for this processor (default value = 1) 16974b9ad928SBarry Smith - mat - the matrices 16984b9ad928SBarry Smith 16994b9ad928SBarry Smith Level: advanced 17004b9ad928SBarry Smith 1701cf739d55SBarry Smith Notes: Call after PCSetUp() (or KSPSetUp()) but before PCApply() (or KSPApply()) and before PCSetUpOnBlocks()) 1702cf739d55SBarry Smith 1703cf739d55SBarry Smith Usually one would use PCSetModifySubmatrices() to change the submatrices in building the preconditioner. 1704cf739d55SBarry Smith 17054b9ad928SBarry Smith .keywords: PC, ASM, set, local, subdomains, additive Schwarz, block Jacobi 17064b9ad928SBarry Smith 17074b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(), 1708cf739d55SBarry Smith PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubdomains(), PCSetModifySubmatrices() 17094b9ad928SBarry Smith @*/ 17107087cfbeSBarry Smith PetscErrorCode PCASMGetLocalSubmatrices(PC pc,PetscInt *n,Mat *mat[]) 17114b9ad928SBarry Smith { 17124b9ad928SBarry Smith PC_ASM *osm; 171392bb6962SLisandro Dalcin PetscErrorCode ierr; 1714ace3abfcSBarry Smith PetscBool match; 17154b9ad928SBarry Smith 17164b9ad928SBarry Smith PetscFunctionBegin; 17170700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 171892bb6962SLisandro Dalcin PetscValidIntPointer(n,2); 171992bb6962SLisandro Dalcin if (mat) PetscValidPointer(mat,3); 1720ce94432eSBarry Smith if (!pc->setupcalled) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call after KSPSetUP() or PCSetUp()."); 1721251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)pc,PCASM,&match);CHKERRQ(ierr); 172292bb6962SLisandro Dalcin if (!match) { 172392bb6962SLisandro Dalcin if (n) *n = 0; 17240298fd71SBarry Smith if (mat) *mat = NULL; 172592bb6962SLisandro Dalcin } else { 17264b9ad928SBarry Smith osm = (PC_ASM*)pc->data; 17274b9ad928SBarry Smith if (n) *n = osm->n_local_true; 17284b9ad928SBarry Smith if (mat) *mat = osm->pmat; 172992bb6962SLisandro Dalcin } 17304b9ad928SBarry Smith PetscFunctionReturn(0); 17314b9ad928SBarry Smith } 1732d709ea83SDmitry Karpeev 1733d709ea83SDmitry Karpeev /*@ 1734d709ea83SDmitry Karpeev PCASMSetDMSubdomains - Indicates whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible. 1735f1ee410cSBarry Smith 1736d709ea83SDmitry Karpeev Logically Collective 1737d709ea83SDmitry Karpeev 1738d709ea83SDmitry Karpeev Input Parameter: 1739d709ea83SDmitry Karpeev + pc - the preconditioner 1740d709ea83SDmitry Karpeev - flg - boolean indicating whether to use subdomains defined by the DM 1741d709ea83SDmitry Karpeev 1742d709ea83SDmitry Karpeev Options Database Key: 1743d709ea83SDmitry Karpeev . -pc_asm_dm_subdomains 1744d709ea83SDmitry Karpeev 1745d709ea83SDmitry Karpeev Level: intermediate 1746d709ea83SDmitry Karpeev 1747d709ea83SDmitry Karpeev Notes: 1748d709ea83SDmitry Karpeev PCASMSetTotalSubdomains() and PCASMSetOverlap() take precedence over PCASMSetDMSubdomains(), 1749d709ea83SDmitry Karpeev so setting either of the first two effectively turns the latter off. 1750d709ea83SDmitry Karpeev 1751d709ea83SDmitry Karpeev .keywords: PC, ASM, DM, set, subdomains, additive Schwarz 1752d709ea83SDmitry Karpeev 1753d709ea83SDmitry Karpeev .seealso: PCASMGetDMSubdomains(), PCASMSetTotalSubdomains(), PCASMSetOverlap() 1754d709ea83SDmitry Karpeev PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubdomains() 1755d709ea83SDmitry Karpeev @*/ 1756d709ea83SDmitry Karpeev PetscErrorCode PCASMSetDMSubdomains(PC pc,PetscBool flg) 1757d709ea83SDmitry Karpeev { 1758d709ea83SDmitry Karpeev PC_ASM *osm = (PC_ASM*)pc->data; 1759d709ea83SDmitry Karpeev PetscErrorCode ierr; 1760d709ea83SDmitry Karpeev PetscBool match; 1761d709ea83SDmitry Karpeev 1762d709ea83SDmitry Karpeev PetscFunctionBegin; 1763d709ea83SDmitry Karpeev PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1764d709ea83SDmitry Karpeev PetscValidLogicalCollectiveBool(pc,flg,2); 1765d709ea83SDmitry Karpeev if (pc->setupcalled) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Not for a setup PC."); 1766d709ea83SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)pc,PCASM,&match);CHKERRQ(ierr); 1767d709ea83SDmitry Karpeev if (match) { 1768d709ea83SDmitry Karpeev osm->dm_subdomains = flg; 1769d709ea83SDmitry Karpeev } 1770d709ea83SDmitry Karpeev PetscFunctionReturn(0); 1771d709ea83SDmitry Karpeev } 1772d709ea83SDmitry Karpeev 1773d709ea83SDmitry Karpeev /*@ 1774d709ea83SDmitry Karpeev PCASMGetDMSubdomains - Returns flag indicating whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible. 1775d709ea83SDmitry Karpeev Not Collective 1776d709ea83SDmitry Karpeev 1777d709ea83SDmitry Karpeev Input Parameter: 1778d709ea83SDmitry Karpeev . pc - the preconditioner 1779d709ea83SDmitry Karpeev 1780d709ea83SDmitry Karpeev Output Parameter: 1781d709ea83SDmitry Karpeev . flg - boolean indicating whether to use subdomains defined by the DM 1782d709ea83SDmitry Karpeev 1783d709ea83SDmitry Karpeev Level: intermediate 1784d709ea83SDmitry Karpeev 1785d709ea83SDmitry Karpeev .keywords: PC, ASM, DM, set, subdomains, additive Schwarz 1786d709ea83SDmitry Karpeev 1787d709ea83SDmitry Karpeev .seealso: PCASMSetDMSubdomains(), PCASMSetTotalSubdomains(), PCASMSetOverlap() 1788d709ea83SDmitry Karpeev PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubdomains() 1789d709ea83SDmitry Karpeev @*/ 1790d709ea83SDmitry Karpeev PetscErrorCode PCASMGetDMSubdomains(PC pc,PetscBool* flg) 1791d709ea83SDmitry Karpeev { 1792d709ea83SDmitry Karpeev PC_ASM *osm = (PC_ASM*)pc->data; 1793d709ea83SDmitry Karpeev PetscErrorCode ierr; 1794d709ea83SDmitry Karpeev PetscBool match; 1795d709ea83SDmitry Karpeev 1796d709ea83SDmitry Karpeev PetscFunctionBegin; 1797d709ea83SDmitry Karpeev PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1798d709ea83SDmitry Karpeev PetscValidPointer(flg,2); 1799d709ea83SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)pc,PCASM,&match);CHKERRQ(ierr); 1800d709ea83SDmitry Karpeev if (match) { 1801d709ea83SDmitry Karpeev if (flg) *flg = osm->dm_subdomains; 1802d709ea83SDmitry Karpeev } 1803d709ea83SDmitry Karpeev PetscFunctionReturn(0); 1804d709ea83SDmitry Karpeev } 180580ec0b7dSPatrick Sanan 180680ec0b7dSPatrick Sanan /*@ 180780ec0b7dSPatrick Sanan PCASMGetSubMatType - Gets the matrix type used for ASM subsolves, as a string. 180880ec0b7dSPatrick Sanan 180980ec0b7dSPatrick Sanan Not Collective 181080ec0b7dSPatrick Sanan 181180ec0b7dSPatrick Sanan Input Parameter: 181280ec0b7dSPatrick Sanan . pc - the PC 181380ec0b7dSPatrick Sanan 181480ec0b7dSPatrick Sanan Output Parameter: 1815f1ee410cSBarry Smith . -pc_asm_sub_mat_type - name of matrix type 181680ec0b7dSPatrick Sanan 181780ec0b7dSPatrick Sanan Level: advanced 181880ec0b7dSPatrick Sanan 181980ec0b7dSPatrick Sanan .keywords: PC, PCASM, MatType, set 182080ec0b7dSPatrick Sanan 182180ec0b7dSPatrick Sanan .seealso: PCASMSetSubMatType(), PCASM, PCSetType(), VecSetType(), MatType, Mat 182280ec0b7dSPatrick Sanan @*/ 182380ec0b7dSPatrick Sanan PetscErrorCode PCASMGetSubMatType(PC pc,MatType *sub_mat_type){ 182480ec0b7dSPatrick Sanan PetscErrorCode ierr; 182580ec0b7dSPatrick Sanan 182680ec0b7dSPatrick Sanan ierr = PetscTryMethod(pc,"PCASMGetSubMatType_C",(PC,MatType*),(pc,sub_mat_type));CHKERRQ(ierr); 182780ec0b7dSPatrick Sanan PetscFunctionReturn(0); 182880ec0b7dSPatrick Sanan } 182980ec0b7dSPatrick Sanan 183080ec0b7dSPatrick Sanan /*@ 183180ec0b7dSPatrick Sanan PCASMSetSubMatType - Set the type of matrix used for ASM subsolves 183280ec0b7dSPatrick Sanan 183380ec0b7dSPatrick Sanan Collective on Mat 183480ec0b7dSPatrick Sanan 183580ec0b7dSPatrick Sanan Input Parameters: 183680ec0b7dSPatrick Sanan + pc - the PC object 183780ec0b7dSPatrick Sanan - sub_mat_type - matrix type 183880ec0b7dSPatrick Sanan 183980ec0b7dSPatrick Sanan Options Database Key: 184080ec0b7dSPatrick 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. 184180ec0b7dSPatrick Sanan 184280ec0b7dSPatrick Sanan Notes: 184380ec0b7dSPatrick Sanan See "${PETSC_DIR}/include/petscmat.h" for available types 184480ec0b7dSPatrick Sanan 184580ec0b7dSPatrick Sanan Level: advanced 184680ec0b7dSPatrick Sanan 184780ec0b7dSPatrick Sanan .keywords: PC, PCASM, MatType, set 184880ec0b7dSPatrick Sanan 184980ec0b7dSPatrick Sanan .seealso: PCASMGetSubMatType(), PCASM, PCSetType(), VecSetType(), MatType, Mat 185080ec0b7dSPatrick Sanan @*/ 185180ec0b7dSPatrick Sanan PetscErrorCode PCASMSetSubMatType(PC pc,MatType sub_mat_type) 185280ec0b7dSPatrick Sanan { 185380ec0b7dSPatrick Sanan PetscErrorCode ierr; 185480ec0b7dSPatrick Sanan 185580ec0b7dSPatrick Sanan ierr = PetscTryMethod(pc,"PCASMSetSubMatType_C",(PC,MatType),(pc,sub_mat_type));CHKERRQ(ierr); 185680ec0b7dSPatrick Sanan PetscFunctionReturn(0); 185780ec0b7dSPatrick Sanan } 1858