14b9ad928SBarry Smith /* 24b9ad928SBarry Smith This file defines an additive Schwarz preconditioner for any Mat implementation. 34b9ad928SBarry Smith 44b9ad928SBarry Smith Note that each processor may have any number of subdomains. But in order to 54b9ad928SBarry Smith deal easily with the VecScatter(), we treat each processor as if it has the 64b9ad928SBarry Smith same number of subdomains. 74b9ad928SBarry Smith 84b9ad928SBarry Smith n - total number of true subdomains on all processors 94b9ad928SBarry Smith n_local_true - actual number of subdomains on this processor 104b9ad928SBarry Smith n_local = maximum over all processors of n_local_true 114b9ad928SBarry Smith */ 124b9ad928SBarry Smith 13910cf402Sprj- #include <../src/ksp/pc/impls/asm/asm.h> /*I "petscpc.h" I*/ 144b9ad928SBarry Smith 156849ba73SBarry Smith static PetscErrorCode PCView_ASM(PC pc,PetscViewer viewer) 164b9ad928SBarry Smith { 1792bb6962SLisandro Dalcin PC_ASM *osm = (PC_ASM*)pc->data; 186849ba73SBarry Smith PetscErrorCode ierr; 1913f74950SBarry Smith PetscMPIInt rank; 204d219a6aSLisandro Dalcin PetscInt i,bsz; 21ace3abfcSBarry Smith PetscBool iascii,isstring; 224b9ad928SBarry Smith PetscViewer sviewer; 234b9ad928SBarry Smith 244b9ad928SBarry Smith PetscFunctionBegin; 25251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 26251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr); 2732077d6dSBarry Smith if (iascii) { 283d03bb48SJed Brown char overlaps[256] = "user-defined overlap",blocks[256] = "total subdomain blocks not yet set"; 298caf3d72SBarry Smith if (osm->overlap >= 0) {ierr = PetscSNPrintf(overlaps,sizeof(overlaps),"amount of overlap = %D",osm->overlap);CHKERRQ(ierr);} 308caf3d72SBarry Smith if (osm->n > 0) {ierr = PetscSNPrintf(blocks,sizeof(blocks),"total subdomain blocks = %D",osm->n);CHKERRQ(ierr);} 31efd4aadfSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," %s, %s\n",blocks,overlaps);CHKERRQ(ierr); 32efd4aadfSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," restriction/interpolation type - %s\n",PCASMTypes[osm->type]);CHKERRQ(ierr); 33e64f0791SPatrick Sanan if (osm->dm_subdomains) {ierr = PetscViewerASCIIPrintf(viewer," Additive Schwarz: using DM to define subdomains\n");CHKERRQ(ierr);} 3412cd4985SMatthew G. Knepley if (osm->loctype != PC_COMPOSITE_ADDITIVE) {ierr = PetscViewerASCIIPrintf(viewer," Additive Schwarz: local solve composition type - %s\n",PCCompositeTypes[osm->loctype]);CHKERRQ(ierr);} 35ce94432eSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);CHKERRQ(ierr); 3692bb6962SLisandro Dalcin if (osm->same_local_solves) { 374d219a6aSLisandro Dalcin if (osm->ksp) { 388d76b567SPatrick Sanan ierr = PetscViewerASCIIPrintf(viewer," Local solver is the same for all blocks, as in the following KSP and PC objects on rank 0:\n");CHKERRQ(ierr); 39c5e4d11fSDmitry Karpeev ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 404d219a6aSLisandro Dalcin if (!rank) { 414b9ad928SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 4292bb6962SLisandro Dalcin ierr = KSPView(osm->ksp[0],sviewer);CHKERRQ(ierr); 434b9ad928SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 444b9ad928SBarry Smith } 45c5e4d11fSDmitry Karpeev ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 464d219a6aSLisandro Dalcin } 474b9ad928SBarry Smith } else { 48c5e4d11fSDmitry Karpeev ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 494d219a6aSLisandro Dalcin ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d] number of local blocks = %D\n",(int)rank,osm->n_local_true);CHKERRQ(ierr); 504d219a6aSLisandro Dalcin ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 514b9ad928SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Local solve info for each block is in the following KSP and PC objects:\n");CHKERRQ(ierr); 524b9ad928SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 534d219a6aSLisandro Dalcin ierr = PetscViewerASCIIPrintf(viewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr); 54c5e4d11fSDmitry Karpeev ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 55dd2fa690SBarry Smith for (i=0; i<osm->n_local_true; i++) { 564d219a6aSLisandro Dalcin ierr = ISGetLocalSize(osm->is[i],&bsz);CHKERRQ(ierr); 574d219a6aSLisandro Dalcin ierr = PetscViewerASCIISynchronizedPrintf(sviewer,"[%d] local block number %D, size = %D\n",(int)rank,i,bsz);CHKERRQ(ierr); 5892bb6962SLisandro Dalcin ierr = KSPView(osm->ksp[i],sviewer);CHKERRQ(ierr); 594d219a6aSLisandro Dalcin ierr = PetscViewerASCIISynchronizedPrintf(sviewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr); 604b9ad928SBarry Smith } 61c5e4d11fSDmitry Karpeev ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 624b9ad928SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 634b9ad928SBarry Smith ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 64c5e4d11fSDmitry Karpeev ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 654b9ad928SBarry Smith } 664b9ad928SBarry Smith } else if (isstring) { 674d219a6aSLisandro Dalcin ierr = PetscViewerStringSPrintf(viewer," blocks=%D, overlap=%D, type=%s",osm->n,osm->overlap,PCASMTypes[osm->type]);CHKERRQ(ierr); 68c5e4d11fSDmitry Karpeev ierr = PetscViewerGetSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 6992bb6962SLisandro Dalcin if (osm->ksp) {ierr = KSPView(osm->ksp[0],sviewer);CHKERRQ(ierr);} 70c5e4d11fSDmitry Karpeev ierr = PetscViewerRestoreSubViewer(viewer,PETSC_COMM_SELF,&sviewer);CHKERRQ(ierr); 714b9ad928SBarry Smith } 724b9ad928SBarry Smith PetscFunctionReturn(0); 734b9ad928SBarry Smith } 744b9ad928SBarry Smith 7592bb6962SLisandro Dalcin static PetscErrorCode PCASMPrintSubdomains(PC pc) 7692bb6962SLisandro Dalcin { 7792bb6962SLisandro Dalcin PC_ASM *osm = (PC_ASM*)pc->data; 7892bb6962SLisandro Dalcin const char *prefix; 7992bb6962SLisandro Dalcin char fname[PETSC_MAX_PATH_LEN+1]; 80643535aeSDmitry Karpeev PetscViewer viewer, sviewer; 81643535aeSDmitry Karpeev char *s; 8292bb6962SLisandro Dalcin PetscInt i,j,nidx; 8392bb6962SLisandro Dalcin const PetscInt *idx; 84643535aeSDmitry Karpeev PetscMPIInt rank, size; 8592bb6962SLisandro Dalcin PetscErrorCode ierr; 86846783a0SBarry Smith 8792bb6962SLisandro Dalcin PetscFunctionBegin; 88ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size);CHKERRQ(ierr); 89ce94432eSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank);CHKERRQ(ierr); 9092bb6962SLisandro Dalcin ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 91589a23caSBarry Smith ierr = PetscOptionsGetString(NULL,prefix,"-pc_asm_print_subdomains",fname,sizeof(fname),NULL);CHKERRQ(ierr); 9292bb6962SLisandro Dalcin if (fname[0] == 0) { ierr = PetscStrcpy(fname,"stdout");CHKERRQ(ierr); }; 93ce94432eSBarry Smith ierr = PetscViewerASCIIOpen(PetscObjectComm((PetscObject)pc),fname,&viewer);CHKERRQ(ierr); 94643535aeSDmitry Karpeev for (i=0; i<osm->n_local; i++) { 95643535aeSDmitry Karpeev if (i < osm->n_local_true) { 9692bb6962SLisandro Dalcin ierr = ISGetLocalSize(osm->is[i],&nidx);CHKERRQ(ierr); 9792bb6962SLisandro Dalcin ierr = ISGetIndices(osm->is[i],&idx);CHKERRQ(ierr); 98643535aeSDmitry Karpeev /* Print to a string viewer; no more than 15 characters per index plus 512 char for the header.*/ 9936a9e3b9SBarry Smith #define len 16*(nidx+1)+512 10036a9e3b9SBarry Smith ierr = PetscMalloc1(len,&s);CHKERRQ(ierr); 10136a9e3b9SBarry Smith ierr = PetscViewerStringOpen(PETSC_COMM_SELF, s, len, &sviewer);CHKERRQ(ierr); 10236a9e3b9SBarry Smith #undef len 103643535aeSDmitry Karpeev ierr = PetscViewerStringSPrintf(sviewer, "[%D:%D] Subdomain %D with overlap:\n", rank, size, i);CHKERRQ(ierr); 10492bb6962SLisandro Dalcin for (j=0; j<nidx; j++) { 105643535aeSDmitry Karpeev ierr = PetscViewerStringSPrintf(sviewer,"%D ",idx[j]);CHKERRQ(ierr); 10692bb6962SLisandro Dalcin } 10792bb6962SLisandro Dalcin ierr = ISRestoreIndices(osm->is[i],&idx);CHKERRQ(ierr); 108643535aeSDmitry Karpeev ierr = PetscViewerStringSPrintf(sviewer,"\n");CHKERRQ(ierr); 109643535aeSDmitry Karpeev ierr = PetscViewerDestroy(&sviewer);CHKERRQ(ierr); 110c5e4d11fSDmitry Karpeev ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 111643535aeSDmitry Karpeev ierr = PetscViewerASCIISynchronizedPrintf(viewer, s);CHKERRQ(ierr); 112643535aeSDmitry Karpeev ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 113c5e4d11fSDmitry Karpeev ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 114643535aeSDmitry Karpeev ierr = PetscFree(s);CHKERRQ(ierr); 1152b691e39SMatthew Knepley if (osm->is_local) { 116643535aeSDmitry Karpeev /* Print to a string viewer; no more than 15 characters per index plus 512 char for the header.*/ 11736a9e3b9SBarry Smith #define len 16*(nidx+1)+512 11836a9e3b9SBarry Smith ierr = PetscMalloc1(len, &s);CHKERRQ(ierr); 11936a9e3b9SBarry Smith ierr = PetscViewerStringOpen(PETSC_COMM_SELF, s, len, &sviewer);CHKERRQ(ierr); 12036a9e3b9SBarry Smith #undef len 12109d011a0SDmitry Karpeev ierr = PetscViewerStringSPrintf(sviewer, "[%D:%D] Subdomain %D without overlap:\n", rank, size, i);CHKERRQ(ierr); 1222b691e39SMatthew Knepley ierr = ISGetLocalSize(osm->is_local[i],&nidx);CHKERRQ(ierr); 1232b691e39SMatthew Knepley ierr = ISGetIndices(osm->is_local[i],&idx);CHKERRQ(ierr); 1242b691e39SMatthew Knepley for (j=0; j<nidx; j++) { 12509d011a0SDmitry Karpeev ierr = PetscViewerStringSPrintf(sviewer,"%D ",idx[j]);CHKERRQ(ierr); 1262b691e39SMatthew Knepley } 1272b691e39SMatthew Knepley ierr = ISRestoreIndices(osm->is_local[i],&idx);CHKERRQ(ierr); 12809d011a0SDmitry Karpeev ierr = PetscViewerStringSPrintf(sviewer,"\n");CHKERRQ(ierr); 129643535aeSDmitry Karpeev ierr = PetscViewerDestroy(&sviewer);CHKERRQ(ierr); 130c5e4d11fSDmitry Karpeev ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 131643535aeSDmitry Karpeev ierr = PetscViewerASCIISynchronizedPrintf(viewer, s);CHKERRQ(ierr); 132643535aeSDmitry Karpeev ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 133c5e4d11fSDmitry Karpeev ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 134643535aeSDmitry Karpeev ierr = PetscFree(s);CHKERRQ(ierr); 135643535aeSDmitry Karpeev } 1362fa5cd67SKarl Rupp } else { 137643535aeSDmitry Karpeev /* Participate in collective viewer calls. */ 138c5e4d11fSDmitry Karpeev ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 139643535aeSDmitry Karpeev ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 140c5e4d11fSDmitry Karpeev ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 141643535aeSDmitry Karpeev /* Assume either all ranks have is_local or none do. */ 142643535aeSDmitry Karpeev if (osm->is_local) { 143c5e4d11fSDmitry Karpeev ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 144643535aeSDmitry Karpeev ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 145c5e4d11fSDmitry Karpeev ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 146643535aeSDmitry Karpeev } 147fdc48646SDmitry Karpeev } 14892bb6962SLisandro Dalcin } 14992bb6962SLisandro Dalcin ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 150fcfd50ebSBarry Smith ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 15192bb6962SLisandro Dalcin PetscFunctionReturn(0); 15292bb6962SLisandro Dalcin } 15392bb6962SLisandro Dalcin 1546849ba73SBarry Smith static PetscErrorCode PCSetUp_ASM(PC pc) 1554b9ad928SBarry Smith { 1564b9ad928SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 1576849ba73SBarry Smith PetscErrorCode ierr; 15857501b6eSBarry Smith PetscBool flg; 15987e86712SBarry Smith PetscInt i,m,m_local; 1604b9ad928SBarry Smith MatReuse scall = MAT_REUSE_MATRIX; 1614b9ad928SBarry Smith IS isl; 1624b9ad928SBarry Smith KSP ksp; 1634b9ad928SBarry Smith PC subpc; 1642dcb1b2aSMatthew Knepley const char *prefix,*pprefix; 16523ce1328SBarry Smith Vec vec; 1660298fd71SBarry Smith DM *domain_dm = NULL; 1674b9ad928SBarry Smith 1684b9ad928SBarry Smith PetscFunctionBegin; 1694b9ad928SBarry Smith if (!pc->setupcalled) { 170265a2120SBarry Smith PetscInt m; 17192bb6962SLisandro Dalcin 1722b837212SDmitry Karpeev /* Note: if subdomains have been set either via PCASMSetTotalSubdomains() or via PCASMSetLocalSubdomains(), osm->n_local_true will not be PETSC_DECIDE */ 1732b837212SDmitry Karpeev if (osm->n_local_true == PETSC_DECIDE) { 17469ca1f37SDmitry Karpeev /* no subdomains given */ 17565db9045SDmitry Karpeev /* try pc->dm first, if allowed */ 176d709ea83SDmitry Karpeev if (osm->dm_subdomains && pc->dm) { 177feb221f8SDmitry Karpeev PetscInt num_domains, d; 178feb221f8SDmitry Karpeev char **domain_names; 1798d4ac253SDmitry Karpeev IS *inner_domain_is, *outer_domain_is; 1808d4ac253SDmitry Karpeev ierr = DMCreateDomainDecomposition(pc->dm, &num_domains, &domain_names, &inner_domain_is, &outer_domain_is, &domain_dm);CHKERRQ(ierr); 181704f0395SPatrick Sanan osm->overlap = -1; /* We do not want to increase the overlap of the IS. 182704f0395SPatrick Sanan A future improvement of this code might allow one to use 183704f0395SPatrick Sanan DM-defined subdomains and also increase the overlap, 184704f0395SPatrick Sanan but that is not currently supported */ 18569ca1f37SDmitry Karpeev if (num_domains) { 1868d4ac253SDmitry Karpeev ierr = PCASMSetLocalSubdomains(pc, num_domains, outer_domain_is, inner_domain_is);CHKERRQ(ierr); 18769ca1f37SDmitry Karpeev } 188feb221f8SDmitry Karpeev for (d = 0; d < num_domains; ++d) { 189a35b7b57SDmitry Karpeev if (domain_names) {ierr = PetscFree(domain_names[d]);CHKERRQ(ierr);} 190a35b7b57SDmitry Karpeev if (inner_domain_is) {ierr = ISDestroy(&inner_domain_is[d]);CHKERRQ(ierr);} 191a35b7b57SDmitry Karpeev if (outer_domain_is) {ierr = ISDestroy(&outer_domain_is[d]);CHKERRQ(ierr);} 192feb221f8SDmitry Karpeev } 193feb221f8SDmitry Karpeev ierr = PetscFree(domain_names);CHKERRQ(ierr); 1948d4ac253SDmitry Karpeev ierr = PetscFree(inner_domain_is);CHKERRQ(ierr); 1958d4ac253SDmitry Karpeev ierr = PetscFree(outer_domain_is);CHKERRQ(ierr); 196feb221f8SDmitry Karpeev } 1972b837212SDmitry Karpeev if (osm->n_local_true == PETSC_DECIDE) { 19869ca1f37SDmitry Karpeev /* still no subdomains; use one subdomain per processor */ 1992b837212SDmitry Karpeev osm->n_local_true = 1; 200feb221f8SDmitry Karpeev } 2012b837212SDmitry Karpeev } 2022b837212SDmitry Karpeev { /* determine the global and max number of subdomains */ 2036ac3741eSJed Brown struct {PetscInt max,sum;} inwork,outwork; 204c10200c1SHong Zhang PetscMPIInt size; 205c10200c1SHong Zhang 2066ac3741eSJed Brown inwork.max = osm->n_local_true; 2076ac3741eSJed Brown inwork.sum = osm->n_local_true; 208367daffbSBarry Smith ierr = MPIU_Allreduce(&inwork,&outwork,1,MPIU_2INT,MPIU_MAXSUM_OP,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr); 2096ac3741eSJed Brown osm->n_local = outwork.max; 2106ac3741eSJed Brown osm->n = outwork.sum; 211c10200c1SHong Zhang 212c10200c1SHong Zhang ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr); 213c10200c1SHong Zhang if (outwork.max == 1 && outwork.sum == size) { 2147dae84e0SHong Zhang /* osm->n_local_true = 1 on all processes, set this option may enable use of optimized MatCreateSubMatrices() implementation */ 215c10200c1SHong Zhang ierr = MatSetOption(pc->pmat,MAT_SUBMAT_SINGLEIS,PETSC_TRUE);CHKERRQ(ierr); 216c10200c1SHong Zhang } 2174b9ad928SBarry Smith } 21878904715SLisandro Dalcin if (!osm->is) { /* create the index sets */ 21978904715SLisandro Dalcin ierr = PCASMCreateSubdomains(pc->pmat,osm->n_local_true,&osm->is);CHKERRQ(ierr); 2204b9ad928SBarry Smith } 221f5234e35SJed Brown if (osm->n_local_true > 1 && !osm->is_local) { 222785e854fSJed Brown ierr = PetscMalloc1(osm->n_local_true,&osm->is_local);CHKERRQ(ierr); 223f5234e35SJed Brown for (i=0; i<osm->n_local_true; i++) { 224f5234e35SJed Brown if (osm->overlap > 0) { /* With positive overlap, osm->is[i] will be modified */ 225f5234e35SJed Brown ierr = ISDuplicate(osm->is[i],&osm->is_local[i]);CHKERRQ(ierr); 226f5234e35SJed Brown ierr = ISCopy(osm->is[i],osm->is_local[i]);CHKERRQ(ierr); 227f5234e35SJed Brown } else { 228f5234e35SJed Brown ierr = PetscObjectReference((PetscObject)osm->is[i]);CHKERRQ(ierr); 229f5234e35SJed Brown osm->is_local[i] = osm->is[i]; 230f5234e35SJed Brown } 231f5234e35SJed Brown } 232f5234e35SJed Brown } 23392bb6962SLisandro Dalcin ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 23490d69ab7SBarry Smith flg = PETSC_FALSE; 235c5929fdfSBarry Smith ierr = PetscOptionsGetBool(NULL,prefix,"-pc_asm_print_subdomains",&flg,NULL);CHKERRQ(ierr); 23692bb6962SLisandro Dalcin if (flg) { ierr = PCASMPrintSubdomains(pc);CHKERRQ(ierr); } 2374b9ad928SBarry Smith 2383d03bb48SJed Brown if (osm->overlap > 0) { 2394b9ad928SBarry Smith /* Extend the "overlapping" regions by a number of steps */ 24092bb6962SLisandro Dalcin ierr = MatIncreaseOverlap(pc->pmat,osm->n_local_true,osm->is,osm->overlap);CHKERRQ(ierr); 2413d03bb48SJed Brown } 2426ed231c7SMatthew Knepley if (osm->sort_indices) { 24392bb6962SLisandro Dalcin for (i=0; i<osm->n_local_true; i++) { 2444b9ad928SBarry Smith ierr = ISSort(osm->is[i]);CHKERRQ(ierr); 2452b691e39SMatthew Knepley if (osm->is_local) { 2462b691e39SMatthew Knepley ierr = ISSort(osm->is_local[i]);CHKERRQ(ierr); 2472b691e39SMatthew Knepley } 2484b9ad928SBarry Smith } 2496ed231c7SMatthew Knepley } 2504b9ad928SBarry Smith 251f6991133SBarry Smith if (!osm->ksp) { 25278904715SLisandro Dalcin /* Create the local solvers */ 253785e854fSJed Brown ierr = PetscMalloc1(osm->n_local_true,&osm->ksp);CHKERRQ(ierr); 254feb221f8SDmitry Karpeev if (domain_dm) { 255feb221f8SDmitry Karpeev ierr = PetscInfo(pc,"Setting up ASM subproblems using the embedded DM\n");CHKERRQ(ierr); 256feb221f8SDmitry Karpeev } 25792bb6962SLisandro Dalcin for (i=0; i<osm->n_local_true; i++) { 2584b9ad928SBarry Smith ierr = KSPCreate(PETSC_COMM_SELF,&ksp);CHKERRQ(ierr); 259422a814eSBarry Smith ierr = KSPSetErrorIfNotConverged(ksp,pc->erroriffailure);CHKERRQ(ierr); 2603bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ksp);CHKERRQ(ierr); 26192bb6962SLisandro Dalcin ierr = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,1);CHKERRQ(ierr); 2624b9ad928SBarry Smith ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr); 2634b9ad928SBarry Smith ierr = KSPGetPC(ksp,&subpc);CHKERRQ(ierr); 2644b9ad928SBarry Smith ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 2654b9ad928SBarry Smith ierr = KSPSetOptionsPrefix(ksp,prefix);CHKERRQ(ierr); 2664b9ad928SBarry Smith ierr = KSPAppendOptionsPrefix(ksp,"sub_");CHKERRQ(ierr); 267feb221f8SDmitry Karpeev if (domain_dm) { 268feb221f8SDmitry Karpeev ierr = KSPSetDM(ksp, domain_dm[i]);CHKERRQ(ierr); 269feb221f8SDmitry Karpeev ierr = KSPSetDMActive(ksp, PETSC_FALSE);CHKERRQ(ierr); 270feb221f8SDmitry Karpeev ierr = DMDestroy(&domain_dm[i]);CHKERRQ(ierr); 271feb221f8SDmitry Karpeev } 2724b9ad928SBarry Smith osm->ksp[i] = ksp; 2734b9ad928SBarry Smith } 274feb221f8SDmitry Karpeev if (domain_dm) { 275feb221f8SDmitry Karpeev ierr = PetscFree(domain_dm);CHKERRQ(ierr); 276feb221f8SDmitry Karpeev } 277f6991133SBarry Smith } 2781dd8081eSeaulisa 279fb745f2cSMatthew G. Knepley ierr = ISConcatenate(PETSC_COMM_SELF, osm->n_local_true, osm->is, &osm->lis);CHKERRQ(ierr); 280fb745f2cSMatthew G. Knepley ierr = ISSortRemoveDups(osm->lis);CHKERRQ(ierr); 281fb745f2cSMatthew G. Knepley ierr = ISGetLocalSize(osm->lis, &m);CHKERRQ(ierr); 282fb745f2cSMatthew G. Knepley ierr = VecCreateSeq(PETSC_COMM_SELF, m, &osm->lx);CHKERRQ(ierr); 283fb745f2cSMatthew G. Knepley ierr = VecDuplicate(osm->lx, &osm->ly);CHKERRQ(ierr); 2841dd8081eSeaulisa 2854b9ad928SBarry Smith scall = MAT_INITIAL_MATRIX; 2864b9ad928SBarry Smith } else { 2874b9ad928SBarry Smith /* 2884b9ad928SBarry Smith Destroy the blocks from the previous iteration 2894b9ad928SBarry Smith */ 290e09e08ccSBarry Smith if (pc->flag == DIFFERENT_NONZERO_PATTERN) { 2914b9ad928SBarry Smith ierr = MatDestroyMatrices(osm->n_local_true,&osm->pmat);CHKERRQ(ierr); 2924b9ad928SBarry Smith scall = MAT_INITIAL_MATRIX; 2934b9ad928SBarry Smith } 2944b9ad928SBarry Smith } 2954b9ad928SBarry Smith 29692bb6962SLisandro Dalcin /* 29792bb6962SLisandro Dalcin Extract out the submatrices 29892bb6962SLisandro Dalcin */ 2997dae84e0SHong Zhang ierr = MatCreateSubMatrices(pc->pmat,osm->n_local_true,osm->is,osm->is,scall,&osm->pmat);CHKERRQ(ierr); 30092bb6962SLisandro Dalcin if (scall == MAT_INITIAL_MATRIX) { 30192bb6962SLisandro Dalcin ierr = PetscObjectGetOptionsPrefix((PetscObject)pc->pmat,&pprefix);CHKERRQ(ierr); 30292bb6962SLisandro Dalcin for (i=0; i<osm->n_local_true; i++) { 3033bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)osm->pmat[i]);CHKERRQ(ierr); 30478904715SLisandro Dalcin ierr = PetscObjectSetOptionsPrefix((PetscObject)osm->pmat[i],pprefix);CHKERRQ(ierr); 30592bb6962SLisandro Dalcin } 30692bb6962SLisandro Dalcin } 30780ec0b7dSPatrick Sanan 30880ec0b7dSPatrick Sanan /* Convert the types of the submatrices (if needbe) */ 30980ec0b7dSPatrick Sanan if (osm->sub_mat_type) { 31080ec0b7dSPatrick Sanan for (i=0; i<osm->n_local_true; i++) { 31180ec0b7dSPatrick Sanan ierr = MatConvert(osm->pmat[i],osm->sub_mat_type,MAT_INPLACE_MATRIX,&(osm->pmat[i]));CHKERRQ(ierr); 31280ec0b7dSPatrick Sanan } 31380ec0b7dSPatrick Sanan } 31480ec0b7dSPatrick Sanan 31580ec0b7dSPatrick Sanan if (!pc->setupcalled) { 31680ec0b7dSPatrick Sanan /* Create the local work vectors (from the local matrices) and scatter contexts */ 3170a545947SLisandro Dalcin ierr = MatCreateVecs(pc->pmat,&vec,NULL);CHKERRQ(ierr); 3185b3c0d42Seaulisa 3191dd8081eSeaulisa 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()"); 3201dd8081eSeaulisa if (osm->is_local && osm->type == PC_ASM_RESTRICT && osm->loctype == PC_COMPOSITE_ADDITIVE) { 3211dd8081eSeaulisa ierr = PetscMalloc1(osm->n_local_true,&osm->lprolongation);CHKERRQ(ierr); 3221dd8081eSeaulisa } 3231dd8081eSeaulisa ierr = PetscMalloc1(osm->n_local_true,&osm->lrestriction);CHKERRQ(ierr); 3241dd8081eSeaulisa ierr = PetscMalloc1(osm->n_local_true,&osm->x);CHKERRQ(ierr); 3251dd8081eSeaulisa ierr = PetscMalloc1(osm->n_local_true,&osm->y);CHKERRQ(ierr); 326347574c9Seaulisa 327347574c9Seaulisa ierr = ISGetLocalSize(osm->lis,&m);CHKERRQ(ierr); 328347574c9Seaulisa ierr = ISCreateStride(PETSC_COMM_SELF,m,0,1,&isl);CHKERRQ(ierr); 3299448b7f1SJunchao Zhang ierr = VecScatterCreate(vec,osm->lis,osm->lx,isl,&osm->restriction);CHKERRQ(ierr); 330347574c9Seaulisa ierr = ISDestroy(&isl);CHKERRQ(ierr); 331347574c9Seaulisa 332347574c9Seaulisa 33380ec0b7dSPatrick Sanan for (i=0; i<osm->n_local_true; ++i) { 3345b3c0d42Seaulisa ISLocalToGlobalMapping ltog; 3355b3c0d42Seaulisa IS isll; 3365b3c0d42Seaulisa const PetscInt *idx_is; 3375b3c0d42Seaulisa PetscInt *idx_lis,nout; 3385b3c0d42Seaulisa 33955817e79SBarry Smith ierr = ISGetLocalSize(osm->is[i],&m);CHKERRQ(ierr); 34055817e79SBarry Smith ierr = MatCreateVecs(osm->pmat[i],&osm->x[i],NULL);CHKERRQ(ierr); 34155817e79SBarry Smith ierr = VecDuplicate(osm->x[i],&osm->y[i]);CHKERRQ(ierr); 34255817e79SBarry Smith 343b0de9f23SBarry Smith /* generate a scatter from ly to y[i] picking all the overlapping is[i] entries */ 3445b3c0d42Seaulisa ierr = ISLocalToGlobalMappingCreateIS(osm->lis,<og);CHKERRQ(ierr); 3455b3c0d42Seaulisa ierr = ISGetLocalSize(osm->is[i],&m);CHKERRQ(ierr); 3465b3c0d42Seaulisa ierr = ISGetIndices(osm->is[i], &idx_is);CHKERRQ(ierr); 3475b3c0d42Seaulisa ierr = PetscMalloc1(m,&idx_lis);CHKERRQ(ierr); 3485b3c0d42Seaulisa ierr = ISGlobalToLocalMappingApply(ltog,IS_GTOLM_DROP,m,idx_is,&nout,idx_lis);CHKERRQ(ierr); 3495b3c0d42Seaulisa if (nout != m) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"is not a subset of lis"); 3505b3c0d42Seaulisa ierr = ISRestoreIndices(osm->is[i], &idx_is);CHKERRQ(ierr); 3515b3c0d42Seaulisa ierr = ISCreateGeneral(PETSC_COMM_SELF,m,idx_lis,PETSC_OWN_POINTER,&isll);CHKERRQ(ierr); 352d8b3b5e3Seaulisa ierr = ISLocalToGlobalMappingDestroy(<og);CHKERRQ(ierr); 3535b3c0d42Seaulisa ierr = ISCreateStride(PETSC_COMM_SELF,m,0,1,&isl);CHKERRQ(ierr); 3549448b7f1SJunchao Zhang ierr = VecScatterCreate(osm->ly,isll,osm->y[i],isl,&osm->lrestriction[i]);CHKERRQ(ierr); 3555b3c0d42Seaulisa ierr = ISDestroy(&isll);CHKERRQ(ierr); 3565b3c0d42Seaulisa ierr = ISDestroy(&isl);CHKERRQ(ierr); 357910cf402Sprj- if (osm->lprolongation) { /* generate a scatter from y[i] to ly picking only the the non-overlapping is_local[i] entries */ 358d8b3b5e3Seaulisa ISLocalToGlobalMapping ltog; 359d8b3b5e3Seaulisa IS isll,isll_local; 360d8b3b5e3Seaulisa const PetscInt *idx_local; 361d8b3b5e3Seaulisa PetscInt *idx1, *idx2, nout; 362d8b3b5e3Seaulisa 363d8b3b5e3Seaulisa ierr = ISGetLocalSize(osm->is_local[i],&m_local);CHKERRQ(ierr); 364d8b3b5e3Seaulisa ierr = ISGetIndices(osm->is_local[i], &idx_local);CHKERRQ(ierr); 365d8b3b5e3Seaulisa 366d8b3b5e3Seaulisa ierr = ISLocalToGlobalMappingCreateIS(osm->is[i],<og);CHKERRQ(ierr); 367d8b3b5e3Seaulisa ierr = PetscMalloc1(m_local,&idx1);CHKERRQ(ierr); 368d8b3b5e3Seaulisa ierr = ISGlobalToLocalMappingApply(ltog,IS_GTOLM_DROP,m_local,idx_local,&nout,idx1);CHKERRQ(ierr); 369d8b3b5e3Seaulisa ierr = ISLocalToGlobalMappingDestroy(<og);CHKERRQ(ierr); 370d8b3b5e3Seaulisa if (nout != m_local) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"is_local not a subset of is"); 371d8b3b5e3Seaulisa ierr = ISCreateGeneral(PETSC_COMM_SELF,m_local,idx1,PETSC_OWN_POINTER,&isll);CHKERRQ(ierr); 372d8b3b5e3Seaulisa 373d8b3b5e3Seaulisa ierr = ISLocalToGlobalMappingCreateIS(osm->lis,<og);CHKERRQ(ierr); 374d8b3b5e3Seaulisa ierr = PetscMalloc1(m_local,&idx2);CHKERRQ(ierr); 375d8b3b5e3Seaulisa ierr = ISGlobalToLocalMappingApply(ltog,IS_GTOLM_DROP,m_local,idx_local,&nout,idx2);CHKERRQ(ierr); 376d8b3b5e3Seaulisa ierr = ISLocalToGlobalMappingDestroy(<og);CHKERRQ(ierr); 377d8b3b5e3Seaulisa if (nout != m_local) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"is_local not a subset of lis"); 378d8b3b5e3Seaulisa ierr = ISCreateGeneral(PETSC_COMM_SELF,m_local,idx2,PETSC_OWN_POINTER,&isll_local);CHKERRQ(ierr); 379d8b3b5e3Seaulisa 380d8b3b5e3Seaulisa ierr = ISRestoreIndices(osm->is_local[i], &idx_local);CHKERRQ(ierr); 3819448b7f1SJunchao Zhang ierr = VecScatterCreate(osm->y[i],isll,osm->ly,isll_local,&osm->lprolongation[i]);CHKERRQ(ierr); 382d8b3b5e3Seaulisa 383d8b3b5e3Seaulisa ierr = ISDestroy(&isll);CHKERRQ(ierr); 384d8b3b5e3Seaulisa ierr = ISDestroy(&isll_local);CHKERRQ(ierr); 385d8b3b5e3Seaulisa } 38680ec0b7dSPatrick Sanan } 38780ec0b7dSPatrick Sanan ierr = VecDestroy(&vec);CHKERRQ(ierr); 38880ec0b7dSPatrick Sanan } 38980ec0b7dSPatrick Sanan 390fb745f2cSMatthew G. Knepley if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE) { 391235cc4ceSMatthew G. Knepley IS *cis; 392235cc4ceSMatthew G. Knepley PetscInt c; 393235cc4ceSMatthew G. Knepley 394235cc4ceSMatthew G. Knepley ierr = PetscMalloc1(osm->n_local_true, &cis);CHKERRQ(ierr); 395235cc4ceSMatthew G. Knepley for (c = 0; c < osm->n_local_true; ++c) cis[c] = osm->lis; 3967dae84e0SHong Zhang ierr = MatCreateSubMatrices(pc->pmat, osm->n_local_true, osm->is, cis, scall, &osm->lmats);CHKERRQ(ierr); 397235cc4ceSMatthew G. Knepley ierr = PetscFree(cis);CHKERRQ(ierr); 398fb745f2cSMatthew G. Knepley } 3994b9ad928SBarry Smith 4004b9ad928SBarry Smith /* Return control to the user so that the submatrices can be modified (e.g., to apply 4014b9ad928SBarry Smith different boundary conditions for the submatrices than for the global problem) */ 40292bb6962SLisandro Dalcin ierr = PCModifySubMatrices(pc,osm->n_local_true,osm->is,osm->is,osm->pmat,pc->modifysubmatricesP);CHKERRQ(ierr); 4034b9ad928SBarry Smith 40492bb6962SLisandro Dalcin /* 40592bb6962SLisandro Dalcin Loop over subdomains putting them into local ksp 40692bb6962SLisandro Dalcin */ 40792bb6962SLisandro Dalcin for (i=0; i<osm->n_local_true; i++) { 40823ee1639SBarry Smith ierr = KSPSetOperators(osm->ksp[i],osm->pmat[i],osm->pmat[i]);CHKERRQ(ierr); 40992bb6962SLisandro Dalcin if (!pc->setupcalled) { 410bf108f30SBarry Smith ierr = KSPSetFromOptions(osm->ksp[i]);CHKERRQ(ierr); 4114b9ad928SBarry Smith } 41292bb6962SLisandro Dalcin } 4134b9ad928SBarry Smith PetscFunctionReturn(0); 4144b9ad928SBarry Smith } 4154b9ad928SBarry Smith 4166849ba73SBarry Smith static PetscErrorCode PCSetUpOnBlocks_ASM(PC pc) 4174b9ad928SBarry Smith { 4184b9ad928SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 4196849ba73SBarry Smith PetscErrorCode ierr; 42013f74950SBarry Smith PetscInt i; 421539c167fSBarry Smith KSPConvergedReason reason; 4224b9ad928SBarry Smith 4234b9ad928SBarry Smith PetscFunctionBegin; 4244b9ad928SBarry Smith for (i=0; i<osm->n_local_true; i++) { 4254b9ad928SBarry Smith ierr = KSPSetUp(osm->ksp[i]);CHKERRQ(ierr); 426539c167fSBarry Smith ierr = KSPGetConvergedReason(osm->ksp[i],&reason);CHKERRQ(ierr); 427c0decd05SBarry Smith if (reason == KSP_DIVERGED_PC_FAILED) { 428261222b2SHong Zhang pc->failedreason = PC_SUBPC_ERROR; 429e0eafd54SHong Zhang } 4304b9ad928SBarry Smith } 4314b9ad928SBarry Smith PetscFunctionReturn(0); 4324b9ad928SBarry Smith } 4334b9ad928SBarry Smith 4346849ba73SBarry Smith static PetscErrorCode PCApply_ASM(PC pc,Vec x,Vec y) 4354b9ad928SBarry Smith { 4364b9ad928SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 4376849ba73SBarry Smith PetscErrorCode ierr; 4381dd8081eSeaulisa PetscInt i,n_local_true = osm->n_local_true; 4394b9ad928SBarry Smith ScatterMode forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE; 4404b9ad928SBarry Smith 4414b9ad928SBarry Smith PetscFunctionBegin; 4424b9ad928SBarry Smith /* 44348e38a8aSPierre Jolivet support for limiting the restriction or interpolation to only local 4444b9ad928SBarry Smith subdomain values (leaving the other values 0). 4454b9ad928SBarry Smith */ 4464b9ad928SBarry Smith if (!(osm->type & PC_ASM_RESTRICT)) { 4474b9ad928SBarry Smith forward = SCATTER_FORWARD_LOCAL; 4484b9ad928SBarry Smith /* have to zero the work RHS since scatter may leave some slots empty */ 4491dd8081eSeaulisa ierr = VecSet(osm->lx, 0.0);CHKERRQ(ierr); 4504b9ad928SBarry Smith } 451347574c9Seaulisa if (!(osm->type & PC_ASM_INTERPOLATE)) { 452347574c9Seaulisa reverse = SCATTER_REVERSE_LOCAL; 453347574c9Seaulisa } 4544b9ad928SBarry Smith 455347574c9Seaulisa if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE || osm->loctype == PC_COMPOSITE_ADDITIVE) { 456b0de9f23SBarry Smith /* zero the global and the local solutions */ 457910cf402Sprj- ierr = VecSet(y, 0.0);CHKERRQ(ierr); 4585b3c0d42Seaulisa ierr = VecSet(osm->ly, 0.0);CHKERRQ(ierr); 459347574c9Seaulisa 46048e38a8aSPierre Jolivet /* copy the global RHS to local RHS including the ghost nodes */ 4611dd8081eSeaulisa ierr = VecScatterBegin(osm->restriction, x, osm->lx, INSERT_VALUES, forward);CHKERRQ(ierr); 4621dd8081eSeaulisa ierr = VecScatterEnd(osm->restriction, x, osm->lx, INSERT_VALUES, forward);CHKERRQ(ierr); 463347574c9Seaulisa 46448e38a8aSPierre Jolivet /* restrict local RHS to the overlapping 0-block RHS */ 4651dd8081eSeaulisa ierr = VecScatterBegin(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward);CHKERRQ(ierr); 4661dd8081eSeaulisa ierr = VecScatterEnd(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward);CHKERRQ(ierr); 467d8b3b5e3Seaulisa 46812cd4985SMatthew G. Knepley /* do the local solves */ 46912cd4985SMatthew G. Knepley for (i = 0; i < n_local_true; ++i) { 470347574c9Seaulisa 471b0de9f23SBarry Smith /* solve the overlapping i-block */ 472fa2bb9feSLisandro Dalcin ierr = PetscLogEventBegin(PC_ApplyOnBlocks, osm->ksp[i], osm->x[i], osm->y[i],0);CHKERRQ(ierr); 47312cd4985SMatthew G. Knepley ierr = KSPSolve(osm->ksp[i], osm->x[i], osm->y[i]);CHKERRQ(ierr); 474c0decd05SBarry Smith ierr = KSPCheckSolve(osm->ksp[i], pc, osm->y[i]);CHKERRQ(ierr); 475fa2bb9feSLisandro Dalcin ierr = PetscLogEventEnd(PC_ApplyOnBlocks, osm->ksp[i], osm->x[i], osm->y[i], 0);CHKERRQ(ierr); 476d8b3b5e3Seaulisa 477910cf402Sprj- if (osm->lprolongation) { /* interpolate the non-overlapping i-block solution to the local solution (only for restrictive additive) */ 4781dd8081eSeaulisa ierr = VecScatterBegin(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward);CHKERRQ(ierr); 4791dd8081eSeaulisa ierr = VecScatterEnd(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward);CHKERRQ(ierr); 48048e38a8aSPierre Jolivet } else { /* interpolate the overlapping i-block solution to the local solution */ 4811dd8081eSeaulisa ierr = VecScatterBegin(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse);CHKERRQ(ierr); 4821dd8081eSeaulisa ierr = VecScatterEnd(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse);CHKERRQ(ierr); 483d8b3b5e3Seaulisa } 484347574c9Seaulisa 485347574c9Seaulisa if (i < n_local_true-1) { 48648e38a8aSPierre Jolivet /* restrict local RHS to the overlapping (i+1)-block RHS */ 4871dd8081eSeaulisa ierr = VecScatterBegin(osm->lrestriction[i+1], osm->lx, osm->x[i+1], INSERT_VALUES, forward);CHKERRQ(ierr); 4881dd8081eSeaulisa ierr = VecScatterEnd(osm->lrestriction[i+1], osm->lx, osm->x[i+1], INSERT_VALUES, forward);CHKERRQ(ierr); 489347574c9Seaulisa 490347574c9Seaulisa if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE) { 4918966356dSPierre Jolivet /* update the overlapping (i+1)-block RHS using the current local solution */ 492347574c9Seaulisa ierr = MatMult(osm->lmats[i+1], osm->ly, osm->y[i+1]);CHKERRQ(ierr); 493347574c9Seaulisa ierr = VecAXPBY(osm->x[i+1],-1.,1., osm->y[i+1]);CHKERRQ(ierr); 4947c3d802fSMatthew G. Knepley } 49512cd4985SMatthew G. Knepley } 49612cd4985SMatthew G. Knepley } 49748e38a8aSPierre Jolivet /* add the local solution to the global solution including the ghost nodes */ 4981dd8081eSeaulisa ierr = VecScatterBegin(osm->restriction, osm->ly, y, ADD_VALUES, reverse);CHKERRQ(ierr); 4991dd8081eSeaulisa ierr = VecScatterEnd(osm->restriction, osm->ly, y, ADD_VALUES, reverse);CHKERRQ(ierr); 500*89c010cfSBarry Smith } else SETERRQ1(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Invalid local composition type: %s", PCCompositeTypes[osm->loctype]); 5014b9ad928SBarry Smith PetscFunctionReturn(0); 5024b9ad928SBarry Smith } 5034b9ad928SBarry Smith 50448e38a8aSPierre Jolivet static PetscErrorCode PCMatApply_ASM(PC pc,Mat X,Mat Y) 50548e38a8aSPierre Jolivet { 50648e38a8aSPierre Jolivet PC_ASM *osm = (PC_ASM*)pc->data; 50748e38a8aSPierre Jolivet Mat Z,W; 50848e38a8aSPierre Jolivet Vec x; 50948e38a8aSPierre Jolivet PetscInt i,m,N; 51048e38a8aSPierre Jolivet ScatterMode forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE; 51148e38a8aSPierre Jolivet PetscErrorCode ierr; 51248e38a8aSPierre Jolivet 51348e38a8aSPierre Jolivet PetscFunctionBegin; 51448e38a8aSPierre Jolivet if (osm->n_local_true > 1) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Not yet implemented"); 51548e38a8aSPierre Jolivet /* 51648e38a8aSPierre Jolivet support for limiting the restriction or interpolation to only local 51748e38a8aSPierre Jolivet subdomain values (leaving the other values 0). 51848e38a8aSPierre Jolivet */ 51948e38a8aSPierre Jolivet if (!(osm->type & PC_ASM_RESTRICT)) { 52048e38a8aSPierre Jolivet forward = SCATTER_FORWARD_LOCAL; 52148e38a8aSPierre Jolivet /* have to zero the work RHS since scatter may leave some slots empty */ 52248e38a8aSPierre Jolivet ierr = VecSet(osm->lx, 0.0);CHKERRQ(ierr); 52348e38a8aSPierre Jolivet } 52448e38a8aSPierre Jolivet if (!(osm->type & PC_ASM_INTERPOLATE)) { 52548e38a8aSPierre Jolivet reverse = SCATTER_REVERSE_LOCAL; 52648e38a8aSPierre Jolivet } 52748e38a8aSPierre Jolivet ierr = VecGetLocalSize(osm->x[0], &m);CHKERRQ(ierr); 52848e38a8aSPierre Jolivet ierr = MatGetSize(X, NULL, &N);CHKERRQ(ierr); 52948e38a8aSPierre Jolivet ierr = MatCreateSeqDense(PETSC_COMM_SELF, m, N, NULL, &Z);CHKERRQ(ierr); 53048e38a8aSPierre Jolivet if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE || osm->loctype == PC_COMPOSITE_ADDITIVE) { 53148e38a8aSPierre Jolivet /* zero the global and the local solutions */ 53248e38a8aSPierre Jolivet ierr = MatZeroEntries(Y);CHKERRQ(ierr); 53348e38a8aSPierre Jolivet ierr = VecSet(osm->ly, 0.0);CHKERRQ(ierr); 53448e38a8aSPierre Jolivet 53548e38a8aSPierre Jolivet for (i = 0; i < N; ++i) { 53648e38a8aSPierre Jolivet ierr = MatDenseGetColumnVecRead(X, i, &x); 53748e38a8aSPierre Jolivet /* copy the global RHS to local RHS including the ghost nodes */ 53848e38a8aSPierre Jolivet ierr = VecScatterBegin(osm->restriction, x, osm->lx, INSERT_VALUES, forward);CHKERRQ(ierr); 53948e38a8aSPierre Jolivet ierr = VecScatterEnd(osm->restriction, x, osm->lx, INSERT_VALUES, forward);CHKERRQ(ierr); 54048e38a8aSPierre Jolivet ierr = MatDenseRestoreColumnVecRead(X, i, &x); 54148e38a8aSPierre Jolivet 54248e38a8aSPierre Jolivet ierr = MatDenseGetColumnVecWrite(Z, i, &x); 54348e38a8aSPierre Jolivet /* restrict local RHS to the overlapping 0-block RHS */ 54448e38a8aSPierre Jolivet ierr = VecScatterBegin(osm->lrestriction[0], osm->lx, x, INSERT_VALUES, forward);CHKERRQ(ierr); 54548e38a8aSPierre Jolivet ierr = VecScatterEnd(osm->lrestriction[0], osm->lx, x, INSERT_VALUES, forward);CHKERRQ(ierr); 54648e38a8aSPierre Jolivet ierr = MatDenseRestoreColumnVecWrite(Z, i, &x); 54748e38a8aSPierre Jolivet } 54848e38a8aSPierre Jolivet ierr = MatCreateSeqDense(PETSC_COMM_SELF, m, N, NULL, &W);CHKERRQ(ierr); 54948e38a8aSPierre Jolivet /* solve the overlapping 0-block */ 55048e38a8aSPierre Jolivet ierr = PetscLogEventBegin(PC_ApplyOnBlocks, osm->ksp[0], Z, W, 0);CHKERRQ(ierr); 55148e38a8aSPierre Jolivet ierr = KSPMatSolve(osm->ksp[0], Z, W);CHKERRQ(ierr); 55248e38a8aSPierre Jolivet ierr = KSPCheckSolve(osm->ksp[0], pc, NULL);CHKERRQ(ierr); 55348e38a8aSPierre Jolivet ierr = PetscLogEventEnd(PC_ApplyOnBlocks, osm->ksp[0], Z, W,0);CHKERRQ(ierr); 55448e38a8aSPierre Jolivet ierr = MatDestroy(&Z);CHKERRQ(ierr); 55548e38a8aSPierre Jolivet 55648e38a8aSPierre Jolivet for (i = 0; i < N; ++i) { 55748e38a8aSPierre Jolivet ierr = VecSet(osm->ly, 0.0);CHKERRQ(ierr); 55848e38a8aSPierre Jolivet ierr = MatDenseGetColumnVecRead(W, i, &x); 55948e38a8aSPierre Jolivet if (osm->lprolongation) { /* interpolate the non-overlapping 0-block solution to the local solution (only for restrictive additive) */ 56048e38a8aSPierre Jolivet ierr = VecScatterBegin(osm->lprolongation[0], x, osm->ly, ADD_VALUES, forward);CHKERRQ(ierr); 56148e38a8aSPierre Jolivet ierr = VecScatterEnd(osm->lprolongation[0], x, osm->ly, ADD_VALUES, forward);CHKERRQ(ierr); 56248e38a8aSPierre Jolivet } else { /* interpolate the overlapping 0-block solution to the local solution */ 56348e38a8aSPierre Jolivet ierr = VecScatterBegin(osm->lrestriction[0], x, osm->ly, ADD_VALUES, reverse);CHKERRQ(ierr); 56448e38a8aSPierre Jolivet ierr = VecScatterEnd(osm->lrestriction[0], x, osm->ly, ADD_VALUES, reverse);CHKERRQ(ierr); 56548e38a8aSPierre Jolivet } 56648e38a8aSPierre Jolivet ierr = MatDenseRestoreColumnVecRead(W, i, &x); 56748e38a8aSPierre Jolivet 56848e38a8aSPierre Jolivet ierr = MatDenseGetColumnVecWrite(Y, i, &x); 56948e38a8aSPierre Jolivet /* add the local solution to the global solution including the ghost nodes */ 57048e38a8aSPierre Jolivet ierr = VecScatterBegin(osm->restriction, osm->ly, x, ADD_VALUES, reverse);CHKERRQ(ierr); 57148e38a8aSPierre Jolivet ierr = VecScatterEnd(osm->restriction, osm->ly, x, ADD_VALUES, reverse);CHKERRQ(ierr); 57248e38a8aSPierre Jolivet ierr = MatDenseRestoreColumnVecWrite(Y, i, &x); 57348e38a8aSPierre Jolivet } 57448e38a8aSPierre Jolivet ierr = MatDestroy(&W);CHKERRQ(ierr); 575*89c010cfSBarry Smith } else SETERRQ1(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Invalid local composition type: %s", PCCompositeTypes[osm->loctype]); 57648e38a8aSPierre Jolivet PetscFunctionReturn(0); 57748e38a8aSPierre Jolivet } 57848e38a8aSPierre Jolivet 5796849ba73SBarry Smith static PetscErrorCode PCApplyTranspose_ASM(PC pc,Vec x,Vec y) 5804b9ad928SBarry Smith { 5814b9ad928SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 5826849ba73SBarry Smith PetscErrorCode ierr; 5831dd8081eSeaulisa PetscInt i,n_local_true = osm->n_local_true; 5844b9ad928SBarry Smith ScatterMode forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE; 5854b9ad928SBarry Smith 5864b9ad928SBarry Smith PetscFunctionBegin; 5874b9ad928SBarry Smith /* 5884b9ad928SBarry Smith Support for limiting the restriction or interpolation to only local 5894b9ad928SBarry Smith subdomain values (leaving the other values 0). 5904b9ad928SBarry Smith 5914b9ad928SBarry Smith Note: these are reversed from the PCApply_ASM() because we are applying the 5924b9ad928SBarry Smith transpose of the three terms 5934b9ad928SBarry Smith */ 594d8b3b5e3Seaulisa 5954b9ad928SBarry Smith if (!(osm->type & PC_ASM_INTERPOLATE)) { 5964b9ad928SBarry Smith forward = SCATTER_FORWARD_LOCAL; 5974b9ad928SBarry Smith /* have to zero the work RHS since scatter may leave some slots empty */ 5981dd8081eSeaulisa ierr = VecSet(osm->lx, 0.0);CHKERRQ(ierr); 5994b9ad928SBarry Smith } 6002fa5cd67SKarl Rupp if (!(osm->type & PC_ASM_RESTRICT)) reverse = SCATTER_REVERSE_LOCAL; 6014b9ad928SBarry Smith 602b0de9f23SBarry Smith /* zero the global and the local solutions */ 603910cf402Sprj- ierr = VecSet(y, 0.0);CHKERRQ(ierr); 604d8b3b5e3Seaulisa ierr = VecSet(osm->ly, 0.0);CHKERRQ(ierr); 605d8b3b5e3Seaulisa 606b0de9f23SBarry Smith /* Copy the global RHS to local RHS including the ghost nodes */ 6071dd8081eSeaulisa ierr = VecScatterBegin(osm->restriction, x, osm->lx, INSERT_VALUES, forward);CHKERRQ(ierr); 6081dd8081eSeaulisa ierr = VecScatterEnd(osm->restriction, x, osm->lx, INSERT_VALUES, forward);CHKERRQ(ierr); 609d8b3b5e3Seaulisa 610b0de9f23SBarry Smith /* Restrict local RHS to the overlapping 0-block RHS */ 6111dd8081eSeaulisa ierr = VecScatterBegin(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward);CHKERRQ(ierr); 6121dd8081eSeaulisa ierr = VecScatterEnd(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward);CHKERRQ(ierr); 613d8b3b5e3Seaulisa 6144b9ad928SBarry Smith /* do the local solves */ 615d8b3b5e3Seaulisa for (i = 0; i < n_local_true; ++i) { 616d8b3b5e3Seaulisa 617b0de9f23SBarry Smith /* solve the overlapping i-block */ 618fa2bb9feSLisandro Dalcin ierr = PetscLogEventBegin(PC_ApplyOnBlocks,osm->ksp[i],osm->x[i],osm->y[i],0);CHKERRQ(ierr); 619e1bcd54cSBarry Smith ierr = KSPSolveTranspose(osm->ksp[i], osm->x[i], osm->y[i]);CHKERRQ(ierr); 620c0decd05SBarry Smith ierr = KSPCheckSolve(osm->ksp[i],pc,osm->y[i]);CHKERRQ(ierr); 621fa2bb9feSLisandro Dalcin ierr = PetscLogEventEnd(PC_ApplyOnBlocks,osm->ksp[i],osm->x[i],osm->y[i],0);CHKERRQ(ierr); 622d8b3b5e3Seaulisa 623910cf402Sprj- if (osm->lprolongation) { /* interpolate the non-overlapping i-block solution to the local solution */ 6241dd8081eSeaulisa ierr = VecScatterBegin(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward);CHKERRQ(ierr); 6251dd8081eSeaulisa ierr = VecScatterEnd(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward);CHKERRQ(ierr); 626910cf402Sprj- } else { /* interpolate the overlapping i-block solution to the local solution */ 6271dd8081eSeaulisa ierr = VecScatterBegin(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse);CHKERRQ(ierr); 6281dd8081eSeaulisa ierr = VecScatterEnd(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse);CHKERRQ(ierr); 6294b9ad928SBarry Smith } 630d8b3b5e3Seaulisa 631d8b3b5e3Seaulisa if (i < n_local_true-1) { 632b0de9f23SBarry Smith /* Restrict local RHS to the overlapping (i+1)-block RHS */ 6331dd8081eSeaulisa ierr = VecScatterBegin(osm->lrestriction[i+1], osm->lx, osm->x[i+1], INSERT_VALUES, forward);CHKERRQ(ierr); 6341dd8081eSeaulisa ierr = VecScatterEnd(osm->lrestriction[i+1], osm->lx, osm->x[i+1], INSERT_VALUES, forward);CHKERRQ(ierr); 6354b9ad928SBarry Smith } 6364b9ad928SBarry Smith } 637b0de9f23SBarry Smith /* Add the local solution to the global solution including the ghost nodes */ 6381dd8081eSeaulisa ierr = VecScatterBegin(osm->restriction, osm->ly, y, ADD_VALUES, reverse);CHKERRQ(ierr); 6391dd8081eSeaulisa ierr = VecScatterEnd(osm->restriction, osm->ly, y, ADD_VALUES, reverse);CHKERRQ(ierr); 640d8b3b5e3Seaulisa 6414b9ad928SBarry Smith PetscFunctionReturn(0); 642d8b3b5e3Seaulisa 6434b9ad928SBarry Smith } 6444b9ad928SBarry Smith 645e91c6855SBarry Smith static PetscErrorCode PCReset_ASM(PC pc) 6464b9ad928SBarry Smith { 6474b9ad928SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 6486849ba73SBarry Smith PetscErrorCode ierr; 64913f74950SBarry Smith PetscInt i; 6504b9ad928SBarry Smith 6514b9ad928SBarry Smith PetscFunctionBegin; 65292bb6962SLisandro Dalcin if (osm->ksp) { 65392bb6962SLisandro Dalcin for (i=0; i<osm->n_local_true; i++) { 654e91c6855SBarry Smith ierr = KSPReset(osm->ksp[i]);CHKERRQ(ierr); 65592bb6962SLisandro Dalcin } 65692bb6962SLisandro Dalcin } 657e09e08ccSBarry Smith if (osm->pmat) { 65892bb6962SLisandro Dalcin if (osm->n_local_true > 0) { 65930a70a9aSHong Zhang ierr = MatDestroySubMatrices(osm->n_local_true,&osm->pmat);CHKERRQ(ierr); 66092bb6962SLisandro Dalcin } 66192bb6962SLisandro Dalcin } 6621dd8081eSeaulisa if (osm->lrestriction) { 6631dd8081eSeaulisa ierr = VecScatterDestroy(&osm->restriction);CHKERRQ(ierr); 6641dd8081eSeaulisa for (i=0; i<osm->n_local_true; i++) { 6651dd8081eSeaulisa ierr = VecScatterDestroy(&osm->lrestriction[i]);CHKERRQ(ierr); 6661dd8081eSeaulisa if (osm->lprolongation) {ierr = VecScatterDestroy(&osm->lprolongation[i]);CHKERRQ(ierr);} 667fcfd50ebSBarry Smith ierr = VecDestroy(&osm->x[i]);CHKERRQ(ierr); 668fcfd50ebSBarry Smith ierr = VecDestroy(&osm->y[i]);CHKERRQ(ierr); 6694b9ad928SBarry Smith } 6701dd8081eSeaulisa ierr = PetscFree(osm->lrestriction);CHKERRQ(ierr); 6711dd8081eSeaulisa if (osm->lprolongation) {ierr = PetscFree(osm->lprolongation);CHKERRQ(ierr);} 67205b42c5fSBarry Smith ierr = PetscFree(osm->x);CHKERRQ(ierr); 67378904715SLisandro Dalcin ierr = PetscFree(osm->y);CHKERRQ(ierr); 6741dd8081eSeaulisa 67592bb6962SLisandro Dalcin } 6762b691e39SMatthew Knepley ierr = PCASMDestroySubdomains(osm->n_local_true,osm->is,osm->is_local);CHKERRQ(ierr); 677fb745f2cSMatthew G. Knepley ierr = ISDestroy(&osm->lis);CHKERRQ(ierr); 678fb745f2cSMatthew G. Knepley ierr = VecDestroy(&osm->lx);CHKERRQ(ierr); 679fb745f2cSMatthew G. Knepley ierr = VecDestroy(&osm->ly);CHKERRQ(ierr); 680347574c9Seaulisa if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE) { 681347574c9Seaulisa ierr = MatDestroyMatrices(osm->n_local_true, &osm->lmats);CHKERRQ(ierr); 682fb745f2cSMatthew G. Knepley } 6832fa5cd67SKarl Rupp 68480ec0b7dSPatrick Sanan ierr = PetscFree(osm->sub_mat_type);CHKERRQ(ierr); 68580ec0b7dSPatrick Sanan 6860a545947SLisandro Dalcin osm->is = NULL; 6870a545947SLisandro Dalcin osm->is_local = NULL; 688e91c6855SBarry Smith PetscFunctionReturn(0); 689e91c6855SBarry Smith } 690e91c6855SBarry Smith 691e91c6855SBarry Smith static PetscErrorCode PCDestroy_ASM(PC pc) 692e91c6855SBarry Smith { 693e91c6855SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 694e91c6855SBarry Smith PetscErrorCode ierr; 695e91c6855SBarry Smith PetscInt i; 696e91c6855SBarry Smith 697e91c6855SBarry Smith PetscFunctionBegin; 698e91c6855SBarry Smith ierr = PCReset_ASM(pc);CHKERRQ(ierr); 699e91c6855SBarry Smith if (osm->ksp) { 700e91c6855SBarry Smith for (i=0; i<osm->n_local_true; i++) { 701fcfd50ebSBarry Smith ierr = KSPDestroy(&osm->ksp[i]);CHKERRQ(ierr); 702e91c6855SBarry Smith } 703e91c6855SBarry Smith ierr = PetscFree(osm->ksp);CHKERRQ(ierr); 704e91c6855SBarry Smith } 705e91c6855SBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 7064b9ad928SBarry Smith PetscFunctionReturn(0); 7074b9ad928SBarry Smith } 7084b9ad928SBarry Smith 7094416b707SBarry Smith static PetscErrorCode PCSetFromOptions_ASM(PetscOptionItems *PetscOptionsObject,PC pc) 7104b9ad928SBarry Smith { 7114b9ad928SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 7126849ba73SBarry Smith PetscErrorCode ierr; 7139dcbbd2bSBarry Smith PetscInt blocks,ovl; 71457501b6eSBarry Smith PetscBool flg; 71592bb6962SLisandro Dalcin PCASMType asmtype; 71612cd4985SMatthew G. Knepley PCCompositeType loctype; 71780ec0b7dSPatrick Sanan char sub_mat_type[256]; 7184b9ad928SBarry Smith 7194b9ad928SBarry Smith PetscFunctionBegin; 720e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"Additive Schwarz options");CHKERRQ(ierr); 721d709ea83SDmitry Karpeev ierr = PetscOptionsBool("-pc_asm_dm_subdomains","Use DMCreateDomainDecomposition() to define subdomains","PCASMSetDMSubdomains",osm->dm_subdomains,&osm->dm_subdomains,&flg);CHKERRQ(ierr); 72290d69ab7SBarry Smith ierr = PetscOptionsInt("-pc_asm_blocks","Number of subdomains","PCASMSetTotalSubdomains",osm->n,&blocks,&flg);CHKERRQ(ierr); 72365db9045SDmitry Karpeev if (flg) { 724f77a5242SKarl Rupp ierr = PCASMSetTotalSubdomains(pc,blocks,NULL,NULL);CHKERRQ(ierr); 725d709ea83SDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 72665db9045SDmitry Karpeev } 727342c94f9SMatthew G. Knepley ierr = PetscOptionsInt("-pc_asm_local_blocks","Number of local subdomains","PCASMSetLocalSubdomains",osm->n_local_true,&blocks,&flg);CHKERRQ(ierr); 728342c94f9SMatthew G. Knepley if (flg) { 729342c94f9SMatthew G. Knepley ierr = PCASMSetLocalSubdomains(pc,blocks,NULL,NULL);CHKERRQ(ierr); 730342c94f9SMatthew G. Knepley osm->dm_subdomains = PETSC_FALSE; 731342c94f9SMatthew G. Knepley } 73290d69ab7SBarry Smith ierr = PetscOptionsInt("-pc_asm_overlap","Number of grid points overlap","PCASMSetOverlap",osm->overlap,&ovl,&flg);CHKERRQ(ierr); 73365db9045SDmitry Karpeev if (flg) { 73465db9045SDmitry Karpeev ierr = PCASMSetOverlap(pc,ovl);CHKERRQ(ierr); 735d709ea83SDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 73665db9045SDmitry Karpeev } 73790d69ab7SBarry Smith flg = PETSC_FALSE; 73890d69ab7SBarry Smith ierr = PetscOptionsEnum("-pc_asm_type","Type of restriction/extension","PCASMSetType",PCASMTypes,(PetscEnum)osm->type,(PetscEnum*)&asmtype,&flg);CHKERRQ(ierr); 73992bb6962SLisandro Dalcin if (flg) {ierr = PCASMSetType(pc,asmtype);CHKERRQ(ierr); } 74012cd4985SMatthew G. Knepley flg = PETSC_FALSE; 74112cd4985SMatthew G. Knepley ierr = PetscOptionsEnum("-pc_asm_local_type","Type of local solver composition","PCASMSetLocalType",PCCompositeTypes,(PetscEnum)osm->loctype,(PetscEnum*)&loctype,&flg);CHKERRQ(ierr); 74212cd4985SMatthew G. Knepley if (flg) {ierr = PCASMSetLocalType(pc,loctype);CHKERRQ(ierr); } 74380ec0b7dSPatrick Sanan ierr = PetscOptionsFList("-pc_asm_sub_mat_type","Subsolve Matrix Type","PCASMSetSubMatType",MatList,NULL,sub_mat_type,256,&flg);CHKERRQ(ierr); 74480ec0b7dSPatrick Sanan if (flg) { 745459726d8SSatish Balay ierr = PCASMSetSubMatType(pc,sub_mat_type);CHKERRQ(ierr); 74680ec0b7dSPatrick Sanan } 7474b9ad928SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 7484b9ad928SBarry Smith PetscFunctionReturn(0); 7494b9ad928SBarry Smith } 7504b9ad928SBarry Smith 7514b9ad928SBarry Smith /*------------------------------------------------------------------------------------*/ 7524b9ad928SBarry Smith 7531e6b0712SBarry Smith static PetscErrorCode PCASMSetLocalSubdomains_ASM(PC pc,PetscInt n,IS is[],IS is_local[]) 7544b9ad928SBarry Smith { 7554b9ad928SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 75692bb6962SLisandro Dalcin PetscErrorCode ierr; 75792bb6962SLisandro Dalcin PetscInt i; 7584b9ad928SBarry Smith 7594b9ad928SBarry Smith PetscFunctionBegin; 760e32f2f54SBarry Smith if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Each process must have 1 or more blocks, n = %D",n); 761ce94432eSBarry 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()."); 762e7e72b3dSBarry Smith 7634b9ad928SBarry Smith if (!pc->setupcalled) { 76492bb6962SLisandro Dalcin if (is) { 76592bb6962SLisandro Dalcin for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)is[i]);CHKERRQ(ierr);} 76692bb6962SLisandro Dalcin } 767832fc9a5SMatthew Knepley if (is_local) { 768832fc9a5SMatthew Knepley for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)is_local[i]);CHKERRQ(ierr);} 769832fc9a5SMatthew Knepley } 7702b691e39SMatthew Knepley ierr = PCASMDestroySubdomains(osm->n_local_true,osm->is,osm->is_local);CHKERRQ(ierr); 7712fa5cd67SKarl Rupp 7724b9ad928SBarry Smith osm->n_local_true = n; 7730a545947SLisandro Dalcin osm->is = NULL; 7740a545947SLisandro Dalcin osm->is_local = NULL; 77592bb6962SLisandro Dalcin if (is) { 776785e854fSJed Brown ierr = PetscMalloc1(n,&osm->is);CHKERRQ(ierr); 7772fa5cd67SKarl Rupp for (i=0; i<n; i++) osm->is[i] = is[i]; 7783d03bb48SJed Brown /* Flag indicating that the user has set overlapping subdomains so PCASM should not increase their size. */ 7793d03bb48SJed Brown osm->overlap = -1; 78092bb6962SLisandro Dalcin } 7812b691e39SMatthew Knepley if (is_local) { 782785e854fSJed Brown ierr = PetscMalloc1(n,&osm->is_local);CHKERRQ(ierr); 7832fa5cd67SKarl Rupp for (i=0; i<n; i++) osm->is_local[i] = is_local[i]; 784a35b7b57SDmitry Karpeev if (!is) { 785785e854fSJed Brown ierr = PetscMalloc1(osm->n_local_true,&osm->is);CHKERRQ(ierr); 786a35b7b57SDmitry Karpeev for (i=0; i<osm->n_local_true; i++) { 787a35b7b57SDmitry Karpeev if (osm->overlap > 0) { /* With positive overlap, osm->is[i] will be modified */ 788a35b7b57SDmitry Karpeev ierr = ISDuplicate(osm->is_local[i],&osm->is[i]);CHKERRQ(ierr); 789a35b7b57SDmitry Karpeev ierr = ISCopy(osm->is_local[i],osm->is[i]);CHKERRQ(ierr); 790a35b7b57SDmitry Karpeev } else { 791a35b7b57SDmitry Karpeev ierr = PetscObjectReference((PetscObject)osm->is_local[i]);CHKERRQ(ierr); 792a35b7b57SDmitry Karpeev osm->is[i] = osm->is_local[i]; 793a35b7b57SDmitry Karpeev } 794a35b7b57SDmitry Karpeev } 795a35b7b57SDmitry Karpeev } 7962b691e39SMatthew Knepley } 7974b9ad928SBarry Smith } 7984b9ad928SBarry Smith PetscFunctionReturn(0); 7994b9ad928SBarry Smith } 8004b9ad928SBarry Smith 8011e6b0712SBarry Smith static PetscErrorCode PCASMSetTotalSubdomains_ASM(PC pc,PetscInt N,IS *is,IS *is_local) 8024b9ad928SBarry Smith { 8034b9ad928SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 8046849ba73SBarry Smith PetscErrorCode ierr; 80513f74950SBarry Smith PetscMPIInt rank,size; 80678904715SLisandro Dalcin PetscInt n; 8074b9ad928SBarry Smith 8084b9ad928SBarry Smith PetscFunctionBegin; 809ce94432eSBarry Smith if (N < 1) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Number of total blocks must be > 0, N = %D",N); 810ce94432eSBarry 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."); 8114b9ad928SBarry Smith 8124b9ad928SBarry Smith /* 813880770e9SJed Brown Split the subdomains equally among all processors 8144b9ad928SBarry Smith */ 815ce94432eSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);CHKERRQ(ierr); 816ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr); 8174b9ad928SBarry Smith n = N/size + ((N % size) > rank); 818e32f2f54SBarry 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); 819e32f2f54SBarry Smith if (pc->setupcalled && n != osm->n_local_true) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"PCASMSetTotalSubdomains() should be called before PCSetUp()."); 8204b9ad928SBarry Smith if (!pc->setupcalled) { 8212b691e39SMatthew Knepley ierr = PCASMDestroySubdomains(osm->n_local_true,osm->is,osm->is_local);CHKERRQ(ierr); 8222fa5cd67SKarl Rupp 8234b9ad928SBarry Smith osm->n_local_true = n; 8240a545947SLisandro Dalcin osm->is = NULL; 8250a545947SLisandro Dalcin osm->is_local = NULL; 8264b9ad928SBarry Smith } 8274b9ad928SBarry Smith PetscFunctionReturn(0); 8284b9ad928SBarry Smith } 8294b9ad928SBarry Smith 8301e6b0712SBarry Smith static PetscErrorCode PCASMSetOverlap_ASM(PC pc,PetscInt ovl) 8314b9ad928SBarry Smith { 83292bb6962SLisandro Dalcin PC_ASM *osm = (PC_ASM*)pc->data; 8334b9ad928SBarry Smith 8344b9ad928SBarry Smith PetscFunctionBegin; 835ce94432eSBarry Smith if (ovl < 0) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap value requested"); 836ce94432eSBarry Smith if (pc->setupcalled && ovl != osm->overlap) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCASMSetOverlap() should be called before PCSetUp()."); 8372fa5cd67SKarl Rupp if (!pc->setupcalled) osm->overlap = ovl; 8384b9ad928SBarry Smith PetscFunctionReturn(0); 8394b9ad928SBarry Smith } 8404b9ad928SBarry Smith 8411e6b0712SBarry Smith static PetscErrorCode PCASMSetType_ASM(PC pc,PCASMType type) 8424b9ad928SBarry Smith { 84392bb6962SLisandro Dalcin PC_ASM *osm = (PC_ASM*)pc->data; 8444b9ad928SBarry Smith 8454b9ad928SBarry Smith PetscFunctionBegin; 8464b9ad928SBarry Smith osm->type = type; 847bf108f30SBarry Smith osm->type_set = PETSC_TRUE; 8484b9ad928SBarry Smith PetscFunctionReturn(0); 8494b9ad928SBarry Smith } 8504b9ad928SBarry Smith 851c60c7ad4SBarry Smith static PetscErrorCode PCASMGetType_ASM(PC pc,PCASMType *type) 852c60c7ad4SBarry Smith { 853c60c7ad4SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 854c60c7ad4SBarry Smith 855c60c7ad4SBarry Smith PetscFunctionBegin; 856c60c7ad4SBarry Smith *type = osm->type; 857c60c7ad4SBarry Smith PetscFunctionReturn(0); 858c60c7ad4SBarry Smith } 859c60c7ad4SBarry Smith 86012cd4985SMatthew G. Knepley static PetscErrorCode PCASMSetLocalType_ASM(PC pc, PCCompositeType type) 86112cd4985SMatthew G. Knepley { 86212cd4985SMatthew G. Knepley PC_ASM *osm = (PC_ASM *) pc->data; 86312cd4985SMatthew G. Knepley 86412cd4985SMatthew G. Knepley PetscFunctionBegin; 865d0ecd4dfSBarry 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"); 86612cd4985SMatthew G. Knepley osm->loctype = type; 86712cd4985SMatthew G. Knepley PetscFunctionReturn(0); 86812cd4985SMatthew G. Knepley } 86912cd4985SMatthew G. Knepley 87012cd4985SMatthew G. Knepley static PetscErrorCode PCASMGetLocalType_ASM(PC pc, PCCompositeType *type) 87112cd4985SMatthew G. Knepley { 87212cd4985SMatthew G. Knepley PC_ASM *osm = (PC_ASM *) pc->data; 87312cd4985SMatthew G. Knepley 87412cd4985SMatthew G. Knepley PetscFunctionBegin; 87512cd4985SMatthew G. Knepley *type = osm->loctype; 87612cd4985SMatthew G. Knepley PetscFunctionReturn(0); 87712cd4985SMatthew G. Knepley } 87812cd4985SMatthew G. Knepley 8791e6b0712SBarry Smith static PetscErrorCode PCASMSetSortIndices_ASM(PC pc,PetscBool doSort) 8806ed231c7SMatthew Knepley { 8816ed231c7SMatthew Knepley PC_ASM *osm = (PC_ASM*)pc->data; 8826ed231c7SMatthew Knepley 8836ed231c7SMatthew Knepley PetscFunctionBegin; 8846ed231c7SMatthew Knepley osm->sort_indices = doSort; 8856ed231c7SMatthew Knepley PetscFunctionReturn(0); 8866ed231c7SMatthew Knepley } 8876ed231c7SMatthew Knepley 8881e6b0712SBarry Smith static PetscErrorCode PCASMGetSubKSP_ASM(PC pc,PetscInt *n_local,PetscInt *first_local,KSP **ksp) 8894b9ad928SBarry Smith { 89092bb6962SLisandro Dalcin PC_ASM *osm = (PC_ASM*)pc->data; 891dfbe8321SBarry Smith PetscErrorCode ierr; 8924b9ad928SBarry Smith 8934b9ad928SBarry Smith PetscFunctionBegin; 89434a84908Sprj- 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"); 8954b9ad928SBarry Smith 8962fa5cd67SKarl Rupp if (n_local) *n_local = osm->n_local_true; 89792bb6962SLisandro Dalcin if (first_local) { 898ce94432eSBarry Smith ierr = MPI_Scan(&osm->n_local_true,first_local,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr); 89992bb6962SLisandro Dalcin *first_local -= osm->n_local_true; 90092bb6962SLisandro Dalcin } 90192bb6962SLisandro Dalcin if (ksp) { 90292bb6962SLisandro Dalcin /* Assume that local solves are now different; not necessarily 90392bb6962SLisandro Dalcin true though! This flag is used only for PCView_ASM() */ 90492bb6962SLisandro Dalcin *ksp = osm->ksp; 90592bb6962SLisandro Dalcin osm->same_local_solves = PETSC_FALSE; 90692bb6962SLisandro Dalcin } 9074b9ad928SBarry Smith PetscFunctionReturn(0); 9084b9ad928SBarry Smith } 9094b9ad928SBarry Smith 91080ec0b7dSPatrick Sanan static PetscErrorCode PCASMGetSubMatType_ASM(PC pc,MatType *sub_mat_type) 91180ec0b7dSPatrick Sanan { 91280ec0b7dSPatrick Sanan PC_ASM *osm = (PC_ASM*)pc->data; 91380ec0b7dSPatrick Sanan 91480ec0b7dSPatrick Sanan PetscFunctionBegin; 91580ec0b7dSPatrick Sanan PetscValidHeaderSpecific(pc,PC_CLASSID,1); 91680ec0b7dSPatrick Sanan PetscValidPointer(sub_mat_type,2); 91780ec0b7dSPatrick Sanan *sub_mat_type = osm->sub_mat_type; 91880ec0b7dSPatrick Sanan PetscFunctionReturn(0); 91980ec0b7dSPatrick Sanan } 92080ec0b7dSPatrick Sanan 92180ec0b7dSPatrick Sanan static PetscErrorCode PCASMSetSubMatType_ASM(PC pc,MatType sub_mat_type) 92280ec0b7dSPatrick Sanan { 92380ec0b7dSPatrick Sanan PetscErrorCode ierr; 92480ec0b7dSPatrick Sanan PC_ASM *osm = (PC_ASM*)pc->data; 92580ec0b7dSPatrick Sanan 92680ec0b7dSPatrick Sanan PetscFunctionBegin; 92780ec0b7dSPatrick Sanan PetscValidHeaderSpecific(pc,PC_CLASSID,1); 92880ec0b7dSPatrick Sanan ierr = PetscFree(osm->sub_mat_type);CHKERRQ(ierr); 92980ec0b7dSPatrick Sanan ierr = PetscStrallocpy(sub_mat_type,(char**)&osm->sub_mat_type);CHKERRQ(ierr); 93080ec0b7dSPatrick Sanan PetscFunctionReturn(0); 93180ec0b7dSPatrick Sanan } 93280ec0b7dSPatrick Sanan 9334b9ad928SBarry Smith /*@C 9341093a601SBarry Smith PCASMSetLocalSubdomains - Sets the local subdomains (for this processor only) for the additive Schwarz preconditioner. 9354b9ad928SBarry Smith 936d083f849SBarry Smith Collective on pc 9374b9ad928SBarry Smith 9384b9ad928SBarry Smith Input Parameters: 9394b9ad928SBarry Smith + pc - the preconditioner context 9404b9ad928SBarry Smith . n - the number of subdomains for this processor (default value = 1) 9418c03b21aSDmitry Karpeev . is - the index set that defines the subdomains for this processor 9420298fd71SBarry Smith (or NULL for PETSc to determine subdomains) 943f1ee410cSBarry 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 944f1ee410cSBarry Smith (or NULL to not provide these) 9454b9ad928SBarry Smith 946342c94f9SMatthew G. Knepley Options Database Key: 947342c94f9SMatthew G. Knepley To set the total number of subdomain blocks rather than specify the 948342c94f9SMatthew G. Knepley index sets, use the option 949342c94f9SMatthew G. Knepley . -pc_asm_local_blocks <blks> - Sets local blocks 950342c94f9SMatthew G. Knepley 9514b9ad928SBarry Smith Notes: 9521093a601SBarry Smith The IS numbering is in the parallel, global numbering of the vector for both is and is_local 9534b9ad928SBarry Smith 9544b9ad928SBarry Smith By default the ASM preconditioner uses 1 block per processor. 9554b9ad928SBarry Smith 9564b9ad928SBarry Smith Use PCASMSetTotalSubdomains() to set the subdomains for all processors. 9574b9ad928SBarry Smith 958f1ee410cSBarry Smith If is_local is provided and PCASMType is PC_ASM_RESTRICT then the solution only over the is_local region is interpolated 959f1ee410cSBarry Smith back to form the global solution (this is the standard restricted additive Schwarz method) 960f1ee410cSBarry Smith 961f1ee410cSBarry 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 962f1ee410cSBarry Smith no code to handle that case. 963f1ee410cSBarry Smith 9644b9ad928SBarry Smith Level: advanced 9654b9ad928SBarry Smith 9664b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(), 967f1ee410cSBarry Smith PCASMCreateSubdomains2D(), PCASMGetLocalSubdomains(), PCASMType, PCASMSetType() 9684b9ad928SBarry Smith @*/ 9697087cfbeSBarry Smith PetscErrorCode PCASMSetLocalSubdomains(PC pc,PetscInt n,IS is[],IS is_local[]) 9704b9ad928SBarry Smith { 9717bb14e67SBarry Smith PetscErrorCode ierr; 9724b9ad928SBarry Smith 9734b9ad928SBarry Smith PetscFunctionBegin; 9740700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 9757bb14e67SBarry Smith ierr = PetscTryMethod(pc,"PCASMSetLocalSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,n,is,is_local));CHKERRQ(ierr); 9764b9ad928SBarry Smith PetscFunctionReturn(0); 9774b9ad928SBarry Smith } 9784b9ad928SBarry Smith 9794b9ad928SBarry Smith /*@C 980feb221f8SDmitry Karpeev PCASMSetTotalSubdomains - Sets the subdomains for all processors for the 9814b9ad928SBarry Smith additive Schwarz preconditioner. Either all or no processors in the 9824b9ad928SBarry Smith PC communicator must call this routine, with the same index sets. 9834b9ad928SBarry Smith 984d083f849SBarry Smith Collective on pc 9854b9ad928SBarry Smith 9864b9ad928SBarry Smith Input Parameters: 9874b9ad928SBarry Smith + pc - the preconditioner context 988feb221f8SDmitry Karpeev . N - the number of subdomains for all processors 989feb221f8SDmitry Karpeev . is - the index sets that define the subdomains for all processors 990dfaa0c5fSBarry Smith (or NULL to ask PETSc to determine the subdomains) 9912b691e39SMatthew Knepley - is_local - the index sets that define the local part of the subdomains for this processor 992f1ee410cSBarry Smith (or NULL to not provide this information) 9934b9ad928SBarry Smith 9944b9ad928SBarry Smith Options Database Key: 9954b9ad928SBarry Smith To set the total number of subdomain blocks rather than specify the 9964b9ad928SBarry Smith index sets, use the option 9974b9ad928SBarry Smith . -pc_asm_blocks <blks> - Sets total blocks 9984b9ad928SBarry Smith 9994b9ad928SBarry Smith Notes: 1000f1ee410cSBarry Smith Currently you cannot use this to set the actual subdomains with the argument is or is_local. 10014b9ad928SBarry Smith 10024b9ad928SBarry Smith By default the ASM preconditioner uses 1 block per processor. 10034b9ad928SBarry Smith 10044b9ad928SBarry Smith These index sets cannot be destroyed until after completion of the 10054b9ad928SBarry Smith linear solves for which the ASM preconditioner is being used. 10064b9ad928SBarry Smith 10074b9ad928SBarry Smith Use PCASMSetLocalSubdomains() to set local subdomains. 10084b9ad928SBarry Smith 10091093a601SBarry Smith The IS numbering is in the parallel, global numbering of the vector for both is and is_local 10101093a601SBarry Smith 10114b9ad928SBarry Smith Level: advanced 10124b9ad928SBarry Smith 10134b9ad928SBarry Smith .seealso: PCASMSetLocalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(), 10144b9ad928SBarry Smith PCASMCreateSubdomains2D() 10154b9ad928SBarry Smith @*/ 10167087cfbeSBarry Smith PetscErrorCode PCASMSetTotalSubdomains(PC pc,PetscInt N,IS is[],IS is_local[]) 10174b9ad928SBarry Smith { 10187bb14e67SBarry Smith PetscErrorCode ierr; 10194b9ad928SBarry Smith 10204b9ad928SBarry Smith PetscFunctionBegin; 10210700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 10227bb14e67SBarry Smith ierr = PetscTryMethod(pc,"PCASMSetTotalSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,N,is,is_local));CHKERRQ(ierr); 10234b9ad928SBarry Smith PetscFunctionReturn(0); 10244b9ad928SBarry Smith } 10254b9ad928SBarry Smith 10264b9ad928SBarry Smith /*@ 10274b9ad928SBarry Smith PCASMSetOverlap - Sets the overlap between a pair of subdomains for the 10284b9ad928SBarry Smith additive Schwarz preconditioner. Either all or no processors in the 1029f1ee410cSBarry Smith PC communicator must call this routine. 10304b9ad928SBarry Smith 1031d083f849SBarry Smith Logically Collective on pc 10324b9ad928SBarry Smith 10334b9ad928SBarry Smith Input Parameters: 10344b9ad928SBarry Smith + pc - the preconditioner context 10354b9ad928SBarry Smith - ovl - the amount of overlap between subdomains (ovl >= 0, default value = 1) 10364b9ad928SBarry Smith 10374b9ad928SBarry Smith Options Database Key: 10384b9ad928SBarry Smith . -pc_asm_overlap <ovl> - Sets overlap 10394b9ad928SBarry Smith 10404b9ad928SBarry Smith Notes: 10414b9ad928SBarry Smith By default the ASM preconditioner uses 1 block per processor. To use 10424b9ad928SBarry Smith multiple blocks per perocessor, see PCASMSetTotalSubdomains() and 10434b9ad928SBarry Smith PCASMSetLocalSubdomains() (and the option -pc_asm_blocks <blks>). 10444b9ad928SBarry Smith 10454b9ad928SBarry Smith The overlap defaults to 1, so if one desires that no additional 10464b9ad928SBarry Smith overlap be computed beyond what may have been set with a call to 10474b9ad928SBarry Smith PCASMSetTotalSubdomains() or PCASMSetLocalSubdomains(), then ovl 10484b9ad928SBarry Smith must be set to be 0. In particular, if one does not explicitly set 10494b9ad928SBarry Smith the subdomains an application code, then all overlap would be computed 10504b9ad928SBarry Smith internally by PETSc, and using an overlap of 0 would result in an ASM 10514b9ad928SBarry Smith variant that is equivalent to the block Jacobi preconditioner. 10524b9ad928SBarry Smith 1053f1ee410cSBarry Smith The default algorithm used by PETSc to increase overlap is fast, but not scalable, 1054f1ee410cSBarry Smith use the option -mat_increase_overlap_scalable when the problem and number of processes is large. 1055f1ee410cSBarry Smith 10564b9ad928SBarry Smith Note that one can define initial index sets with any overlap via 1057f1ee410cSBarry Smith PCASMSetLocalSubdomains(); the routine 10584b9ad928SBarry Smith PCASMSetOverlap() merely allows PETSc to extend that overlap further 10594b9ad928SBarry Smith if desired. 10604b9ad928SBarry Smith 10614b9ad928SBarry Smith Level: intermediate 10624b9ad928SBarry Smith 10634b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetLocalSubdomains(), PCASMGetSubKSP(), 1064f1ee410cSBarry Smith PCASMCreateSubdomains2D(), PCASMGetLocalSubdomains(), MatIncreaseOverlap() 10654b9ad928SBarry Smith @*/ 10667087cfbeSBarry Smith PetscErrorCode PCASMSetOverlap(PC pc,PetscInt ovl) 10674b9ad928SBarry Smith { 10687bb14e67SBarry Smith PetscErrorCode ierr; 10694b9ad928SBarry Smith 10704b9ad928SBarry Smith PetscFunctionBegin; 10710700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1072c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,ovl,2); 10737bb14e67SBarry Smith ierr = PetscTryMethod(pc,"PCASMSetOverlap_C",(PC,PetscInt),(pc,ovl));CHKERRQ(ierr); 10744b9ad928SBarry Smith PetscFunctionReturn(0); 10754b9ad928SBarry Smith } 10764b9ad928SBarry Smith 10774b9ad928SBarry Smith /*@ 10784b9ad928SBarry Smith PCASMSetType - Sets the type of restriction and interpolation used 10794b9ad928SBarry Smith for local problems in the additive Schwarz method. 10804b9ad928SBarry Smith 1081d083f849SBarry Smith Logically Collective on pc 10824b9ad928SBarry Smith 10834b9ad928SBarry Smith Input Parameters: 10844b9ad928SBarry Smith + pc - the preconditioner context 10854b9ad928SBarry Smith - type - variant of ASM, one of 10864b9ad928SBarry Smith .vb 10874b9ad928SBarry Smith PC_ASM_BASIC - full interpolation and restriction 10884b9ad928SBarry Smith PC_ASM_RESTRICT - full restriction, local processor interpolation 10894b9ad928SBarry Smith PC_ASM_INTERPOLATE - full interpolation, local processor restriction 10904b9ad928SBarry Smith PC_ASM_NONE - local processor restriction and interpolation 10914b9ad928SBarry Smith .ve 10924b9ad928SBarry Smith 10934b9ad928SBarry Smith Options Database Key: 10944b9ad928SBarry Smith . -pc_asm_type [basic,restrict,interpolate,none] - Sets ASM type 10954b9ad928SBarry Smith 109695452b02SPatrick Sanan Notes: 109795452b02SPatrick Sanan if the is_local arguments are passed to PCASMSetLocalSubdomains() then they are used when PC_ASM_RESTRICT has been selected 1098f1ee410cSBarry Smith to limit the local processor interpolation 1099f1ee410cSBarry Smith 11004b9ad928SBarry Smith Level: intermediate 11014b9ad928SBarry Smith 11024b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetTotalSubdomains(), PCASMGetSubKSP(), 1103f1ee410cSBarry Smith PCASMCreateSubdomains2D(), PCASMType, PCASMSetLocalType(), PCASMGetLocalType() 11044b9ad928SBarry Smith @*/ 11057087cfbeSBarry Smith PetscErrorCode PCASMSetType(PC pc,PCASMType type) 11064b9ad928SBarry Smith { 11077bb14e67SBarry Smith PetscErrorCode ierr; 11084b9ad928SBarry Smith 11094b9ad928SBarry Smith PetscFunctionBegin; 11100700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1111c5eb9154SBarry Smith PetscValidLogicalCollectiveEnum(pc,type,2); 11127bb14e67SBarry Smith ierr = PetscTryMethod(pc,"PCASMSetType_C",(PC,PCASMType),(pc,type));CHKERRQ(ierr); 11134b9ad928SBarry Smith PetscFunctionReturn(0); 11144b9ad928SBarry Smith } 11154b9ad928SBarry Smith 1116c60c7ad4SBarry Smith /*@ 1117c60c7ad4SBarry Smith PCASMGetType - Gets the type of restriction and interpolation used 1118c60c7ad4SBarry Smith for local problems in the additive Schwarz method. 1119c60c7ad4SBarry Smith 1120d083f849SBarry Smith Logically Collective on pc 1121c60c7ad4SBarry Smith 1122c60c7ad4SBarry Smith Input Parameter: 1123c60c7ad4SBarry Smith . pc - the preconditioner context 1124c60c7ad4SBarry Smith 1125c60c7ad4SBarry Smith Output Parameter: 1126c60c7ad4SBarry Smith . type - variant of ASM, one of 1127c60c7ad4SBarry Smith 1128c60c7ad4SBarry Smith .vb 1129c60c7ad4SBarry Smith PC_ASM_BASIC - full interpolation and restriction 1130c60c7ad4SBarry Smith PC_ASM_RESTRICT - full restriction, local processor interpolation 1131c60c7ad4SBarry Smith PC_ASM_INTERPOLATE - full interpolation, local processor restriction 1132c60c7ad4SBarry Smith PC_ASM_NONE - local processor restriction and interpolation 1133c60c7ad4SBarry Smith .ve 1134c60c7ad4SBarry Smith 1135c60c7ad4SBarry Smith Options Database Key: 1136c60c7ad4SBarry Smith . -pc_asm_type [basic,restrict,interpolate,none] - Sets ASM type 1137c60c7ad4SBarry Smith 1138c60c7ad4SBarry Smith Level: intermediate 1139c60c7ad4SBarry Smith 1140c60c7ad4SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetTotalSubdomains(), PCASMGetSubKSP(), 1141f1ee410cSBarry Smith PCASMCreateSubdomains2D(), PCASMType, PCASMSetType(), PCASMSetLocalType(), PCASMGetLocalType() 1142c60c7ad4SBarry Smith @*/ 1143c60c7ad4SBarry Smith PetscErrorCode PCASMGetType(PC pc,PCASMType *type) 1144c60c7ad4SBarry Smith { 1145c60c7ad4SBarry Smith PetscErrorCode ierr; 1146c60c7ad4SBarry Smith 1147c60c7ad4SBarry Smith PetscFunctionBegin; 1148c60c7ad4SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1149c60c7ad4SBarry Smith ierr = PetscUseMethod(pc,"PCASMGetType_C",(PC,PCASMType*),(pc,type));CHKERRQ(ierr); 1150c60c7ad4SBarry Smith PetscFunctionReturn(0); 1151c60c7ad4SBarry Smith } 1152c60c7ad4SBarry Smith 115312cd4985SMatthew G. Knepley /*@ 115412cd4985SMatthew G. Knepley PCASMSetLocalType - Sets the type of composition used for local problems in the additive Schwarz method. 115512cd4985SMatthew G. Knepley 1156d083f849SBarry Smith Logically Collective on pc 115712cd4985SMatthew G. Knepley 115812cd4985SMatthew G. Knepley Input Parameters: 115912cd4985SMatthew G. Knepley + pc - the preconditioner context 116012cd4985SMatthew G. Knepley - type - type of composition, one of 116112cd4985SMatthew G. Knepley .vb 116212cd4985SMatthew G. Knepley PC_COMPOSITE_ADDITIVE - local additive combination 116312cd4985SMatthew G. Knepley PC_COMPOSITE_MULTIPLICATIVE - local multiplicative combination 116412cd4985SMatthew G. Knepley .ve 116512cd4985SMatthew G. Knepley 116612cd4985SMatthew G. Knepley Options Database Key: 116712cd4985SMatthew G. Knepley . -pc_asm_local_type [additive,multiplicative] - Sets local solver composition type 116812cd4985SMatthew G. Knepley 116912cd4985SMatthew G. Knepley Level: intermediate 117012cd4985SMatthew G. Knepley 1171f1ee410cSBarry Smith .seealso: PCASMSetType(), PCASMGetType(), PCASMGetLocalType(), PCASM, PCASMType, PCASMSetType(), PCASMGetType(), PCCompositeType 117212cd4985SMatthew G. Knepley @*/ 117312cd4985SMatthew G. Knepley PetscErrorCode PCASMSetLocalType(PC pc, PCCompositeType type) 117412cd4985SMatthew G. Knepley { 117512cd4985SMatthew G. Knepley PetscErrorCode ierr; 117612cd4985SMatthew G. Knepley 117712cd4985SMatthew G. Knepley PetscFunctionBegin; 117812cd4985SMatthew G. Knepley PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 117912cd4985SMatthew G. Knepley PetscValidLogicalCollectiveEnum(pc, type, 2); 118012cd4985SMatthew G. Knepley ierr = PetscTryMethod(pc, "PCASMSetLocalType_C", (PC, PCCompositeType), (pc, type));CHKERRQ(ierr); 118112cd4985SMatthew G. Knepley PetscFunctionReturn(0); 118212cd4985SMatthew G. Knepley } 118312cd4985SMatthew G. Knepley 118412cd4985SMatthew G. Knepley /*@ 118512cd4985SMatthew G. Knepley PCASMGetLocalType - Gets the type of composition used for local problems in the additive Schwarz method. 118612cd4985SMatthew G. Knepley 1187d083f849SBarry Smith Logically Collective on pc 118812cd4985SMatthew G. Knepley 118912cd4985SMatthew G. Knepley Input Parameter: 119012cd4985SMatthew G. Knepley . pc - the preconditioner context 119112cd4985SMatthew G. Knepley 119212cd4985SMatthew G. Knepley Output Parameter: 119312cd4985SMatthew G. Knepley . type - type of composition, one of 119412cd4985SMatthew G. Knepley .vb 119512cd4985SMatthew G. Knepley PC_COMPOSITE_ADDITIVE - local additive combination 119612cd4985SMatthew G. Knepley PC_COMPOSITE_MULTIPLICATIVE - local multiplicative combination 119712cd4985SMatthew G. Knepley .ve 119812cd4985SMatthew G. Knepley 119912cd4985SMatthew G. Knepley Options Database Key: 120012cd4985SMatthew G. Knepley . -pc_asm_local_type [additive,multiplicative] - Sets local solver composition type 120112cd4985SMatthew G. Knepley 120212cd4985SMatthew G. Knepley Level: intermediate 120312cd4985SMatthew G. Knepley 1204f1ee410cSBarry Smith .seealso: PCASMSetType(), PCASMGetType(), PCASMSetLocalType(), PCASMCreate(), PCASMType, PCASMSetType(), PCASMGetType(), PCCompositeType 120512cd4985SMatthew G. Knepley @*/ 120612cd4985SMatthew G. Knepley PetscErrorCode PCASMGetLocalType(PC pc, PCCompositeType *type) 120712cd4985SMatthew G. Knepley { 120812cd4985SMatthew G. Knepley PetscErrorCode ierr; 120912cd4985SMatthew G. Knepley 121012cd4985SMatthew G. Knepley PetscFunctionBegin; 121112cd4985SMatthew G. Knepley PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 121212cd4985SMatthew G. Knepley PetscValidPointer(type, 2); 121312cd4985SMatthew G. Knepley ierr = PetscUseMethod(pc, "PCASMGetLocalType_C", (PC, PCCompositeType *), (pc, type));CHKERRQ(ierr); 121412cd4985SMatthew G. Knepley PetscFunctionReturn(0); 121512cd4985SMatthew G. Knepley } 121612cd4985SMatthew G. Knepley 12176ed231c7SMatthew Knepley /*@ 12186ed231c7SMatthew Knepley PCASMSetSortIndices - Determines whether subdomain indices are sorted. 12196ed231c7SMatthew Knepley 1220d083f849SBarry Smith Logically Collective on pc 12216ed231c7SMatthew Knepley 12226ed231c7SMatthew Knepley Input Parameters: 12236ed231c7SMatthew Knepley + pc - the preconditioner context 12246ed231c7SMatthew Knepley - doSort - sort the subdomain indices 12256ed231c7SMatthew Knepley 12266ed231c7SMatthew Knepley Level: intermediate 12276ed231c7SMatthew Knepley 12286ed231c7SMatthew Knepley .seealso: PCASMSetLocalSubdomains(), PCASMSetTotalSubdomains(), PCASMGetSubKSP(), 12296ed231c7SMatthew Knepley PCASMCreateSubdomains2D() 12306ed231c7SMatthew Knepley @*/ 12317087cfbeSBarry Smith PetscErrorCode PCASMSetSortIndices(PC pc,PetscBool doSort) 12326ed231c7SMatthew Knepley { 12337bb14e67SBarry Smith PetscErrorCode ierr; 12346ed231c7SMatthew Knepley 12356ed231c7SMatthew Knepley PetscFunctionBegin; 12360700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1237acfcf0e5SJed Brown PetscValidLogicalCollectiveBool(pc,doSort,2); 12387bb14e67SBarry Smith ierr = PetscTryMethod(pc,"PCASMSetSortIndices_C",(PC,PetscBool),(pc,doSort));CHKERRQ(ierr); 12396ed231c7SMatthew Knepley PetscFunctionReturn(0); 12406ed231c7SMatthew Knepley } 12416ed231c7SMatthew Knepley 12424b9ad928SBarry Smith /*@C 12434b9ad928SBarry Smith PCASMGetSubKSP - Gets the local KSP contexts for all blocks on 12444b9ad928SBarry Smith this processor. 12454b9ad928SBarry Smith 1246d083f849SBarry Smith Collective on pc iff first_local is requested 12474b9ad928SBarry Smith 12484b9ad928SBarry Smith Input Parameter: 12494b9ad928SBarry Smith . pc - the preconditioner context 12504b9ad928SBarry Smith 12514b9ad928SBarry Smith Output Parameters: 12520298fd71SBarry Smith + n_local - the number of blocks on this processor or NULL 12530298fd71SBarry Smith . first_local - the global number of the first block on this processor or NULL, 12540298fd71SBarry Smith all processors must request or all must pass NULL 12554b9ad928SBarry Smith - ksp - the array of KSP contexts 12564b9ad928SBarry Smith 12574b9ad928SBarry Smith Note: 1258d29017ddSJed Brown After PCASMGetSubKSP() the array of KSPes is not to be freed. 12594b9ad928SBarry Smith 12604b9ad928SBarry Smith You must call KSPSetUp() before calling PCASMGetSubKSP(). 12614b9ad928SBarry Smith 1262d29017ddSJed Brown Fortran note: 12632bf68e3eSBarry 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. 1264d29017ddSJed Brown 12654b9ad928SBarry Smith Level: advanced 12664b9ad928SBarry Smith 12674b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetTotalSubdomains(), PCASMSetOverlap(), 12684b9ad928SBarry Smith PCASMCreateSubdomains2D(), 12694b9ad928SBarry Smith @*/ 12707087cfbeSBarry Smith PetscErrorCode PCASMGetSubKSP(PC pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp[]) 12714b9ad928SBarry Smith { 12727bb14e67SBarry Smith PetscErrorCode ierr; 12734b9ad928SBarry Smith 12744b9ad928SBarry Smith PetscFunctionBegin; 12750700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 12767bb14e67SBarry Smith ierr = PetscUseMethod(pc,"PCASMGetSubKSP_C",(PC,PetscInt*,PetscInt*,KSP **),(pc,n_local,first_local,ksp));CHKERRQ(ierr); 12774b9ad928SBarry Smith PetscFunctionReturn(0); 12784b9ad928SBarry Smith } 12794b9ad928SBarry Smith 12804b9ad928SBarry Smith /* -------------------------------------------------------------------------------------*/ 12814b9ad928SBarry Smith /*MC 12823b09bd56SBarry Smith PCASM - Use the (restricted) additive Schwarz method, each block is (approximately) solved with 12834b9ad928SBarry Smith its own KSP object. 12844b9ad928SBarry Smith 12854b9ad928SBarry Smith Options Database Keys: 128649517cdeSBarry Smith + -pc_asm_blocks <blks> - Sets total blocks 12874b9ad928SBarry Smith . -pc_asm_overlap <ovl> - Sets overlap 1288f1ee410cSBarry Smith . -pc_asm_type [basic,restrict,interpolate,none] - Sets ASM type, default is restrict 1289f1ee410cSBarry Smith - -pc_asm_local_type [additive, multiplicative] - Sets ASM type, default is additive 12904b9ad928SBarry Smith 12913b09bd56SBarry Smith IMPORTANT: If you run with, for example, 3 blocks on 1 processor or 3 blocks on 3 processors you 12923b09bd56SBarry Smith will get a different convergence rate due to the default option of -pc_asm_type restrict. Use 12933b09bd56SBarry Smith -pc_asm_type basic to use the standard ASM. 12943b09bd56SBarry Smith 129595452b02SPatrick Sanan Notes: 129695452b02SPatrick Sanan Each processor can have one or more blocks, but a block cannot be shared by more 1297f1ee410cSBarry Smith than one processor. Use PCGASM for subdomains shared by multiple processes. Defaults to one block per processor. 12984b9ad928SBarry Smith 12993b09bd56SBarry Smith To set options on the solvers for each block append -sub_ to all the KSP, and PC 1300d7ee0231SBarry Smith options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly 13014b9ad928SBarry Smith 1302a8c7a070SBarry Smith To set the options on the solvers separate for each block call PCASMGetSubKSP() 13034b9ad928SBarry Smith and set the options directly on the resulting KSP object (you can access its PC 13044b9ad928SBarry Smith with KSPGetPC()) 13054b9ad928SBarry Smith 13064b9ad928SBarry Smith Level: beginner 13074b9ad928SBarry Smith 1308c582cd25SBarry Smith References: 130996a0c994SBarry Smith + 1. - M Dryja, OB Widlund, An additive variant of the Schwarz alternating method for the case of many subregions 131096a0c994SBarry Smith Courant Institute, New York University Technical report 13116d33885cSprj- - 2. - Barry Smith, Petter Bjorstad, and William Gropp, Domain Decompositions: Parallel Multilevel Methods for Elliptic Partial Differential Equations, 131296a0c994SBarry Smith Cambridge University Press. 1313c582cd25SBarry Smith 13144b9ad928SBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 1315f1ee410cSBarry Smith PCBJACOBI, PCASMGetSubKSP(), PCASMSetLocalSubdomains(), PCASMType, PCASMGetType(), PCASMSetLocalType(), PCASMGetLocalType() 131634a84908Sprj- PCASMSetTotalSubdomains(), PCSetModifySubMatrices(), PCASMSetOverlap(), PCASMSetType(), PCCompositeType 1317e09e08ccSBarry Smith 13184b9ad928SBarry Smith M*/ 13194b9ad928SBarry Smith 13208cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_ASM(PC pc) 13214b9ad928SBarry Smith { 1322dfbe8321SBarry Smith PetscErrorCode ierr; 13234b9ad928SBarry Smith PC_ASM *osm; 13244b9ad928SBarry Smith 13254b9ad928SBarry Smith PetscFunctionBegin; 1326b00a9115SJed Brown ierr = PetscNewLog(pc,&osm);CHKERRQ(ierr); 13272fa5cd67SKarl Rupp 13284b9ad928SBarry Smith osm->n = PETSC_DECIDE; 13294b9ad928SBarry Smith osm->n_local = 0; 13302b837212SDmitry Karpeev osm->n_local_true = PETSC_DECIDE; 13314b9ad928SBarry Smith osm->overlap = 1; 13320a545947SLisandro Dalcin osm->ksp = NULL; 13330a545947SLisandro Dalcin osm->restriction = NULL; 13340a545947SLisandro Dalcin osm->lprolongation = NULL; 13350a545947SLisandro Dalcin osm->lrestriction = NULL; 13360a545947SLisandro Dalcin osm->x = NULL; 13370a545947SLisandro Dalcin osm->y = NULL; 13380a545947SLisandro Dalcin osm->is = NULL; 13390a545947SLisandro Dalcin osm->is_local = NULL; 13400a545947SLisandro Dalcin osm->mat = NULL; 13410a545947SLisandro Dalcin osm->pmat = NULL; 13424b9ad928SBarry Smith osm->type = PC_ASM_RESTRICT; 134312cd4985SMatthew G. Knepley osm->loctype = PC_COMPOSITE_ADDITIVE; 13444b9ad928SBarry Smith osm->same_local_solves = PETSC_TRUE; 13456ed231c7SMatthew Knepley osm->sort_indices = PETSC_TRUE; 1346d709ea83SDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 134780ec0b7dSPatrick Sanan osm->sub_mat_type = NULL; 13484b9ad928SBarry Smith 134992bb6962SLisandro Dalcin pc->data = (void*)osm; 13504b9ad928SBarry Smith pc->ops->apply = PCApply_ASM; 135148e38a8aSPierre Jolivet pc->ops->matapply = PCMatApply_ASM; 13524b9ad928SBarry Smith pc->ops->applytranspose = PCApplyTranspose_ASM; 13534b9ad928SBarry Smith pc->ops->setup = PCSetUp_ASM; 1354e91c6855SBarry Smith pc->ops->reset = PCReset_ASM; 13554b9ad928SBarry Smith pc->ops->destroy = PCDestroy_ASM; 13564b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_ASM; 13574b9ad928SBarry Smith pc->ops->setuponblocks = PCSetUpOnBlocks_ASM; 13584b9ad928SBarry Smith pc->ops->view = PCView_ASM; 13590a545947SLisandro Dalcin pc->ops->applyrichardson = NULL; 13604b9ad928SBarry Smith 1361bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetLocalSubdomains_C",PCASMSetLocalSubdomains_ASM);CHKERRQ(ierr); 1362bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetTotalSubdomains_C",PCASMSetTotalSubdomains_ASM);CHKERRQ(ierr); 1363bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetOverlap_C",PCASMSetOverlap_ASM);CHKERRQ(ierr); 1364bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetType_C",PCASMSetType_ASM);CHKERRQ(ierr); 1365c60c7ad4SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMGetType_C",PCASMGetType_ASM);CHKERRQ(ierr); 136612cd4985SMatthew G. Knepley ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetLocalType_C",PCASMSetLocalType_ASM);CHKERRQ(ierr); 136712cd4985SMatthew G. Knepley ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMGetLocalType_C",PCASMGetLocalType_ASM);CHKERRQ(ierr); 1368bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetSortIndices_C",PCASMSetSortIndices_ASM);CHKERRQ(ierr); 1369bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMGetSubKSP_C",PCASMGetSubKSP_ASM);CHKERRQ(ierr); 137080ec0b7dSPatrick Sanan ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMGetSubMatType_C",PCASMGetSubMatType_ASM);CHKERRQ(ierr); 137180ec0b7dSPatrick Sanan ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetSubMatType_C",PCASMSetSubMatType_ASM);CHKERRQ(ierr); 13724b9ad928SBarry Smith PetscFunctionReturn(0); 13734b9ad928SBarry Smith } 13744b9ad928SBarry Smith 137592bb6962SLisandro Dalcin /*@C 137692bb6962SLisandro Dalcin PCASMCreateSubdomains - Creates the index sets for the overlapping Schwarz 137792bb6962SLisandro Dalcin preconditioner for a any problem on a general grid. 137892bb6962SLisandro Dalcin 137992bb6962SLisandro Dalcin Collective 138092bb6962SLisandro Dalcin 138192bb6962SLisandro Dalcin Input Parameters: 138292bb6962SLisandro Dalcin + A - The global matrix operator 138392bb6962SLisandro Dalcin - n - the number of local blocks 138492bb6962SLisandro Dalcin 138592bb6962SLisandro Dalcin Output Parameters: 138692bb6962SLisandro Dalcin . outis - the array of index sets defining the subdomains 138792bb6962SLisandro Dalcin 138892bb6962SLisandro Dalcin Level: advanced 138992bb6962SLisandro Dalcin 13907d6bfa3bSBarry Smith Note: this generates nonoverlapping subdomains; the PCASM will generate the overlap 13917d6bfa3bSBarry Smith from these if you use PCASMSetLocalSubdomains() 13927d6bfa3bSBarry Smith 13937d6bfa3bSBarry Smith In the Fortran version you must provide the array outis[] already allocated of length n. 13947d6bfa3bSBarry Smith 139592bb6962SLisandro Dalcin .seealso: PCASMSetLocalSubdomains(), PCASMDestroySubdomains() 139692bb6962SLisandro Dalcin @*/ 13977087cfbeSBarry Smith PetscErrorCode PCASMCreateSubdomains(Mat A, PetscInt n, IS* outis[]) 139892bb6962SLisandro Dalcin { 139992bb6962SLisandro Dalcin MatPartitioning mpart; 140092bb6962SLisandro Dalcin const char *prefix; 140192bb6962SLisandro Dalcin PetscInt i,j,rstart,rend,bs; 1402976e8c5aSLisandro Dalcin PetscBool hasop, isbaij = PETSC_FALSE,foundpart = PETSC_FALSE; 14030298fd71SBarry Smith Mat Ad = NULL, adj; 140492bb6962SLisandro Dalcin IS ispart,isnumb,*is; 140592bb6962SLisandro Dalcin PetscErrorCode ierr; 140692bb6962SLisandro Dalcin 140792bb6962SLisandro Dalcin PetscFunctionBegin; 14080700a824SBarry Smith PetscValidHeaderSpecific(A,MAT_CLASSID,1); 140992bb6962SLisandro Dalcin PetscValidPointer(outis,3); 1410c1235816SBarry Smith if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of local blocks must be > 0, n = %D",n); 141192bb6962SLisandro Dalcin 141292bb6962SLisandro Dalcin /* Get prefix, row distribution, and block size */ 141392bb6962SLisandro Dalcin ierr = MatGetOptionsPrefix(A,&prefix);CHKERRQ(ierr); 141492bb6962SLisandro Dalcin ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 141592bb6962SLisandro Dalcin ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr); 141665e19b50SBarry 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); 141765e19b50SBarry Smith 141892bb6962SLisandro Dalcin /* Get diagonal block from matrix if possible */ 1419976e8c5aSLisandro Dalcin ierr = MatHasOperation(A,MATOP_GET_DIAGONAL_BLOCK,&hasop);CHKERRQ(ierr); 1420976e8c5aSLisandro Dalcin if (hasop) { 142111bd1e4dSLisandro Dalcin ierr = MatGetDiagonalBlock(A,&Ad);CHKERRQ(ierr); 142292bb6962SLisandro Dalcin } 142392bb6962SLisandro Dalcin if (Ad) { 1424b9e7e5c1SBarry Smith ierr = PetscObjectBaseTypeCompare((PetscObject)Ad,MATSEQBAIJ,&isbaij);CHKERRQ(ierr); 1425b9e7e5c1SBarry Smith if (!isbaij) {ierr = PetscObjectBaseTypeCompare((PetscObject)Ad,MATSEQSBAIJ,&isbaij);CHKERRQ(ierr);} 142692bb6962SLisandro Dalcin } 142792bb6962SLisandro Dalcin if (Ad && n > 1) { 1428ace3abfcSBarry Smith PetscBool match,done; 142992bb6962SLisandro Dalcin /* Try to setup a good matrix partitioning if available */ 143092bb6962SLisandro Dalcin ierr = MatPartitioningCreate(PETSC_COMM_SELF,&mpart);CHKERRQ(ierr); 143192bb6962SLisandro Dalcin ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr); 143292bb6962SLisandro Dalcin ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr); 1433251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGCURRENT,&match);CHKERRQ(ierr); 143492bb6962SLisandro Dalcin if (!match) { 1435251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGSQUARE,&match);CHKERRQ(ierr); 143692bb6962SLisandro Dalcin } 143792bb6962SLisandro Dalcin if (!match) { /* assume a "good" partitioner is available */ 14381a83f524SJed Brown PetscInt na; 14391a83f524SJed Brown const PetscInt *ia,*ja; 144092bb6962SLisandro Dalcin ierr = MatGetRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr); 144192bb6962SLisandro Dalcin if (done) { 144292bb6962SLisandro Dalcin /* Build adjacency matrix by hand. Unfortunately a call to 144392bb6962SLisandro Dalcin MatConvert(Ad,MATMPIADJ,MAT_INITIAL_MATRIX,&adj) will 144492bb6962SLisandro Dalcin remove the block-aij structure and we cannot expect 144592bb6962SLisandro Dalcin MatPartitioning to split vertices as we need */ 14460a545947SLisandro Dalcin PetscInt i,j,len,nnz,cnt,*iia=NULL,*jja=NULL; 14471a83f524SJed Brown const PetscInt *row; 144892bb6962SLisandro Dalcin nnz = 0; 144992bb6962SLisandro Dalcin for (i=0; i<na; i++) { /* count number of nonzeros */ 145092bb6962SLisandro Dalcin len = ia[i+1] - ia[i]; 145192bb6962SLisandro Dalcin row = ja + ia[i]; 145292bb6962SLisandro Dalcin for (j=0; j<len; j++) { 145392bb6962SLisandro Dalcin if (row[j] == i) { /* don't count diagonal */ 145492bb6962SLisandro Dalcin len--; break; 145592bb6962SLisandro Dalcin } 145692bb6962SLisandro Dalcin } 145792bb6962SLisandro Dalcin nnz += len; 145892bb6962SLisandro Dalcin } 1459854ce69bSBarry Smith ierr = PetscMalloc1(na+1,&iia);CHKERRQ(ierr); 1460854ce69bSBarry Smith ierr = PetscMalloc1(nnz,&jja);CHKERRQ(ierr); 146192bb6962SLisandro Dalcin nnz = 0; 146292bb6962SLisandro Dalcin iia[0] = 0; 146392bb6962SLisandro Dalcin for (i=0; i<na; i++) { /* fill adjacency */ 146492bb6962SLisandro Dalcin cnt = 0; 146592bb6962SLisandro Dalcin len = ia[i+1] - ia[i]; 146692bb6962SLisandro Dalcin row = ja + ia[i]; 146792bb6962SLisandro Dalcin for (j=0; j<len; j++) { 146892bb6962SLisandro Dalcin if (row[j] != i) { /* if not diagonal */ 146992bb6962SLisandro Dalcin jja[nnz+cnt++] = row[j]; 147092bb6962SLisandro Dalcin } 147192bb6962SLisandro Dalcin } 147292bb6962SLisandro Dalcin nnz += cnt; 147392bb6962SLisandro Dalcin iia[i+1] = nnz; 147492bb6962SLisandro Dalcin } 147592bb6962SLisandro Dalcin /* Partitioning of the adjacency matrix */ 14760298fd71SBarry Smith ierr = MatCreateMPIAdj(PETSC_COMM_SELF,na,na,iia,jja,NULL,&adj);CHKERRQ(ierr); 147792bb6962SLisandro Dalcin ierr = MatPartitioningSetAdjacency(mpart,adj);CHKERRQ(ierr); 147892bb6962SLisandro Dalcin ierr = MatPartitioningSetNParts(mpart,n);CHKERRQ(ierr); 147992bb6962SLisandro Dalcin ierr = MatPartitioningApply(mpart,&ispart);CHKERRQ(ierr); 148092bb6962SLisandro Dalcin ierr = ISPartitioningToNumbering(ispart,&isnumb);CHKERRQ(ierr); 1481fcfd50ebSBarry Smith ierr = MatDestroy(&adj);CHKERRQ(ierr); 148292bb6962SLisandro Dalcin foundpart = PETSC_TRUE; 148392bb6962SLisandro Dalcin } 148492bb6962SLisandro Dalcin ierr = MatRestoreRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr); 148592bb6962SLisandro Dalcin } 1486fcfd50ebSBarry Smith ierr = MatPartitioningDestroy(&mpart);CHKERRQ(ierr); 148792bb6962SLisandro Dalcin } 148892bb6962SLisandro Dalcin 1489785e854fSJed Brown ierr = PetscMalloc1(n,&is);CHKERRQ(ierr); 149092bb6962SLisandro Dalcin *outis = is; 149192bb6962SLisandro Dalcin 149292bb6962SLisandro Dalcin if (!foundpart) { 149392bb6962SLisandro Dalcin 149492bb6962SLisandro Dalcin /* Partitioning by contiguous chunks of rows */ 149592bb6962SLisandro Dalcin 149692bb6962SLisandro Dalcin PetscInt mbs = (rend-rstart)/bs; 149792bb6962SLisandro Dalcin PetscInt start = rstart; 149892bb6962SLisandro Dalcin for (i=0; i<n; i++) { 149992bb6962SLisandro Dalcin PetscInt count = (mbs/n + ((mbs % n) > i)) * bs; 150092bb6962SLisandro Dalcin ierr = ISCreateStride(PETSC_COMM_SELF,count,start,1,&is[i]);CHKERRQ(ierr); 150192bb6962SLisandro Dalcin start += count; 150292bb6962SLisandro Dalcin } 150392bb6962SLisandro Dalcin 150492bb6962SLisandro Dalcin } else { 150592bb6962SLisandro Dalcin 150692bb6962SLisandro Dalcin /* Partitioning by adjacency of diagonal block */ 150792bb6962SLisandro Dalcin 150892bb6962SLisandro Dalcin const PetscInt *numbering; 150992bb6962SLisandro Dalcin PetscInt *count,nidx,*indices,*newidx,start=0; 151092bb6962SLisandro Dalcin /* Get node count in each partition */ 1511785e854fSJed Brown ierr = PetscMalloc1(n,&count);CHKERRQ(ierr); 151292bb6962SLisandro Dalcin ierr = ISPartitioningCount(ispart,n,count);CHKERRQ(ierr); 151392bb6962SLisandro Dalcin if (isbaij && bs > 1) { /* adjust for the block-aij case */ 151492bb6962SLisandro Dalcin for (i=0; i<n; i++) count[i] *= bs; 151592bb6962SLisandro Dalcin } 151692bb6962SLisandro Dalcin /* Build indices from node numbering */ 151792bb6962SLisandro Dalcin ierr = ISGetLocalSize(isnumb,&nidx);CHKERRQ(ierr); 1518785e854fSJed Brown ierr = PetscMalloc1(nidx,&indices);CHKERRQ(ierr); 151992bb6962SLisandro Dalcin for (i=0; i<nidx; i++) indices[i] = i; /* needs to be initialized */ 152092bb6962SLisandro Dalcin ierr = ISGetIndices(isnumb,&numbering);CHKERRQ(ierr); 152192bb6962SLisandro Dalcin ierr = PetscSortIntWithPermutation(nidx,numbering,indices);CHKERRQ(ierr); 152292bb6962SLisandro Dalcin ierr = ISRestoreIndices(isnumb,&numbering);CHKERRQ(ierr); 152392bb6962SLisandro Dalcin if (isbaij && bs > 1) { /* adjust for the block-aij case */ 1524785e854fSJed Brown ierr = PetscMalloc1(nidx*bs,&newidx);CHKERRQ(ierr); 15252fa5cd67SKarl Rupp for (i=0; i<nidx; i++) { 15262fa5cd67SKarl Rupp for (j=0; j<bs; j++) newidx[i*bs+j] = indices[i]*bs + j; 15272fa5cd67SKarl Rupp } 152892bb6962SLisandro Dalcin ierr = PetscFree(indices);CHKERRQ(ierr); 152992bb6962SLisandro Dalcin nidx *= bs; 153092bb6962SLisandro Dalcin indices = newidx; 153192bb6962SLisandro Dalcin } 153292bb6962SLisandro Dalcin /* Shift to get global indices */ 153392bb6962SLisandro Dalcin for (i=0; i<nidx; i++) indices[i] += rstart; 153492bb6962SLisandro Dalcin 153592bb6962SLisandro Dalcin /* Build the index sets for each block */ 153692bb6962SLisandro Dalcin for (i=0; i<n; i++) { 153770b3c8c7SBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,count[i],&indices[start],PETSC_COPY_VALUES,&is[i]);CHKERRQ(ierr); 153892bb6962SLisandro Dalcin ierr = ISSort(is[i]);CHKERRQ(ierr); 153992bb6962SLisandro Dalcin start += count[i]; 154092bb6962SLisandro Dalcin } 154192bb6962SLisandro Dalcin 15423bf036e2SBarry Smith ierr = PetscFree(count);CHKERRQ(ierr); 15433bf036e2SBarry Smith ierr = PetscFree(indices);CHKERRQ(ierr); 1544fcfd50ebSBarry Smith ierr = ISDestroy(&isnumb);CHKERRQ(ierr); 1545fcfd50ebSBarry Smith ierr = ISDestroy(&ispart);CHKERRQ(ierr); 154692bb6962SLisandro Dalcin 154792bb6962SLisandro Dalcin } 154892bb6962SLisandro Dalcin PetscFunctionReturn(0); 154992bb6962SLisandro Dalcin } 155092bb6962SLisandro Dalcin 155192bb6962SLisandro Dalcin /*@C 155292bb6962SLisandro Dalcin PCASMDestroySubdomains - Destroys the index sets created with 155392bb6962SLisandro Dalcin PCASMCreateSubdomains(). Should be called after setting subdomains 155492bb6962SLisandro Dalcin with PCASMSetLocalSubdomains(). 155592bb6962SLisandro Dalcin 155692bb6962SLisandro Dalcin Collective 155792bb6962SLisandro Dalcin 155892bb6962SLisandro Dalcin Input Parameters: 155992bb6962SLisandro Dalcin + n - the number of index sets 15602b691e39SMatthew Knepley . is - the array of index sets 15610298fd71SBarry Smith - is_local - the array of local index sets, can be NULL 156292bb6962SLisandro Dalcin 156392bb6962SLisandro Dalcin Level: advanced 156492bb6962SLisandro Dalcin 156592bb6962SLisandro Dalcin .seealso: PCASMCreateSubdomains(), PCASMSetLocalSubdomains() 156692bb6962SLisandro Dalcin @*/ 15677087cfbeSBarry Smith PetscErrorCode PCASMDestroySubdomains(PetscInt n, IS is[], IS is_local[]) 156892bb6962SLisandro Dalcin { 156992bb6962SLisandro Dalcin PetscInt i; 157092bb6962SLisandro Dalcin PetscErrorCode ierr; 15715fd66863SKarl Rupp 157292bb6962SLisandro Dalcin PetscFunctionBegin; 1573a35b7b57SDmitry Karpeev if (n <= 0) PetscFunctionReturn(0); 1574a35b7b57SDmitry Karpeev if (is) { 157592bb6962SLisandro Dalcin PetscValidPointer(is,2); 1576fcfd50ebSBarry Smith for (i=0; i<n; i++) { ierr = ISDestroy(&is[i]);CHKERRQ(ierr); } 157792bb6962SLisandro Dalcin ierr = PetscFree(is);CHKERRQ(ierr); 1578a35b7b57SDmitry Karpeev } 15792b691e39SMatthew Knepley if (is_local) { 15802b691e39SMatthew Knepley PetscValidPointer(is_local,3); 1581fcfd50ebSBarry Smith for (i=0; i<n; i++) { ierr = ISDestroy(&is_local[i]);CHKERRQ(ierr); } 15822b691e39SMatthew Knepley ierr = PetscFree(is_local);CHKERRQ(ierr); 15832b691e39SMatthew Knepley } 158492bb6962SLisandro Dalcin PetscFunctionReturn(0); 158592bb6962SLisandro Dalcin } 158692bb6962SLisandro Dalcin 15874b9ad928SBarry Smith /*@ 15884b9ad928SBarry Smith PCASMCreateSubdomains2D - Creates the index sets for the overlapping Schwarz 15894b9ad928SBarry Smith preconditioner for a two-dimensional problem on a regular grid. 15904b9ad928SBarry Smith 15914b9ad928SBarry Smith Not Collective 15924b9ad928SBarry Smith 15934b9ad928SBarry Smith Input Parameters: 15944b9ad928SBarry Smith + m, n - the number of mesh points in the x and y directions 15954b9ad928SBarry Smith . M, N - the number of subdomains in the x and y directions 15964b9ad928SBarry Smith . dof - degrees of freedom per node 15974b9ad928SBarry Smith - overlap - overlap in mesh lines 15984b9ad928SBarry Smith 15994b9ad928SBarry Smith Output Parameters: 16004b9ad928SBarry Smith + Nsub - the number of subdomains created 16013d03bb48SJed Brown . is - array of index sets defining overlapping (if overlap > 0) subdomains 16023d03bb48SJed Brown - is_local - array of index sets defining non-overlapping subdomains 16034b9ad928SBarry Smith 16044b9ad928SBarry Smith Note: 16054b9ad928SBarry Smith Presently PCAMSCreateSubdomains2d() is valid only for sequential 16064b9ad928SBarry Smith preconditioners. More general related routines are 16074b9ad928SBarry Smith PCASMSetTotalSubdomains() and PCASMSetLocalSubdomains(). 16084b9ad928SBarry Smith 16094b9ad928SBarry Smith Level: advanced 16104b9ad928SBarry Smith 16114b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetLocalSubdomains(), PCASMGetSubKSP(), 16124b9ad928SBarry Smith PCASMSetOverlap() 16134b9ad928SBarry Smith @*/ 16147087cfbeSBarry Smith PetscErrorCode PCASMCreateSubdomains2D(PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt dof,PetscInt overlap,PetscInt *Nsub,IS **is,IS **is_local) 16154b9ad928SBarry Smith { 16163d03bb48SJed Brown PetscInt i,j,height,width,ystart,xstart,yleft,yright,xleft,xright,loc_outer; 16176849ba73SBarry Smith PetscErrorCode ierr; 161813f74950SBarry Smith PetscInt nidx,*idx,loc,ii,jj,count; 16194b9ad928SBarry Smith 16204b9ad928SBarry Smith PetscFunctionBegin; 1621e32f2f54SBarry Smith if (dof != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP," "); 16224b9ad928SBarry Smith 16234b9ad928SBarry Smith *Nsub = N*M; 1624854ce69bSBarry Smith ierr = PetscMalloc1(*Nsub,is);CHKERRQ(ierr); 1625854ce69bSBarry Smith ierr = PetscMalloc1(*Nsub,is_local);CHKERRQ(ierr); 16264b9ad928SBarry Smith ystart = 0; 16273d03bb48SJed Brown loc_outer = 0; 16284b9ad928SBarry Smith for (i=0; i<N; i++) { 16294b9ad928SBarry Smith height = n/N + ((n % N) > i); /* height of subdomain */ 1630e32f2f54SBarry Smith if (height < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many N subdomains for mesh dimension n"); 16314b9ad928SBarry Smith yleft = ystart - overlap; if (yleft < 0) yleft = 0; 16324b9ad928SBarry Smith yright = ystart + height + overlap; if (yright > n) yright = n; 16334b9ad928SBarry Smith xstart = 0; 16344b9ad928SBarry Smith for (j=0; j<M; j++) { 16354b9ad928SBarry Smith width = m/M + ((m % M) > j); /* width of subdomain */ 1636e32f2f54SBarry Smith if (width < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many M subdomains for mesh dimension m"); 16374b9ad928SBarry Smith xleft = xstart - overlap; if (xleft < 0) xleft = 0; 16384b9ad928SBarry Smith xright = xstart + width + overlap; if (xright > m) xright = m; 16394b9ad928SBarry Smith nidx = (xright - xleft)*(yright - yleft); 1640785e854fSJed Brown ierr = PetscMalloc1(nidx,&idx);CHKERRQ(ierr); 16414b9ad928SBarry Smith loc = 0; 16424b9ad928SBarry Smith for (ii=yleft; ii<yright; ii++) { 16434b9ad928SBarry Smith count = m*ii + xleft; 16442fa5cd67SKarl Rupp for (jj=xleft; jj<xright; jj++) idx[loc++] = count++; 16454b9ad928SBarry Smith } 164670b3c8c7SBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,nidx,idx,PETSC_COPY_VALUES,(*is)+loc_outer);CHKERRQ(ierr); 16473d03bb48SJed Brown if (overlap == 0) { 16483d03bb48SJed Brown ierr = PetscObjectReference((PetscObject)(*is)[loc_outer]);CHKERRQ(ierr); 16492fa5cd67SKarl Rupp 16503d03bb48SJed Brown (*is_local)[loc_outer] = (*is)[loc_outer]; 16513d03bb48SJed Brown } else { 16523d03bb48SJed Brown for (loc=0,ii=ystart; ii<ystart+height; ii++) { 16533d03bb48SJed Brown for (jj=xstart; jj<xstart+width; jj++) { 16543d03bb48SJed Brown idx[loc++] = m*ii + jj; 16553d03bb48SJed Brown } 16563d03bb48SJed Brown } 165770b3c8c7SBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,loc,idx,PETSC_COPY_VALUES,*is_local+loc_outer);CHKERRQ(ierr); 16583d03bb48SJed Brown } 16594b9ad928SBarry Smith ierr = PetscFree(idx);CHKERRQ(ierr); 16604b9ad928SBarry Smith xstart += width; 16613d03bb48SJed Brown loc_outer++; 16624b9ad928SBarry Smith } 16634b9ad928SBarry Smith ystart += height; 16644b9ad928SBarry Smith } 16654b9ad928SBarry Smith for (i=0; i<*Nsub; i++) { ierr = ISSort((*is)[i]);CHKERRQ(ierr); } 16664b9ad928SBarry Smith PetscFunctionReturn(0); 16674b9ad928SBarry Smith } 16684b9ad928SBarry Smith 16694b9ad928SBarry Smith /*@C 16704b9ad928SBarry Smith PCASMGetLocalSubdomains - Gets the local subdomains (for this processor 16714b9ad928SBarry Smith only) for the additive Schwarz preconditioner. 16724b9ad928SBarry Smith 1673ad4df100SBarry Smith Not Collective 16744b9ad928SBarry Smith 16754b9ad928SBarry Smith Input Parameter: 16764b9ad928SBarry Smith . pc - the preconditioner context 16774b9ad928SBarry Smith 16784b9ad928SBarry Smith Output Parameters: 16794b9ad928SBarry Smith + n - the number of subdomains for this processor (default value = 1) 16802b691e39SMatthew Knepley . is - the index sets that define the subdomains for this processor 16810298fd71SBarry Smith - is_local - the index sets that define the local part of the subdomains for this processor (can be NULL) 16824b9ad928SBarry Smith 16834b9ad928SBarry Smith 16844b9ad928SBarry Smith Notes: 16854b9ad928SBarry Smith The IS numbering is in the parallel, global numbering of the vector. 16864b9ad928SBarry Smith 16874b9ad928SBarry Smith Level: advanced 16884b9ad928SBarry Smith 16894b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(), 16904b9ad928SBarry Smith PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubmatrices() 16914b9ad928SBarry Smith @*/ 16927087cfbeSBarry Smith PetscErrorCode PCASMGetLocalSubdomains(PC pc,PetscInt *n,IS *is[],IS *is_local[]) 16934b9ad928SBarry Smith { 16942a808120SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 169592bb6962SLisandro Dalcin PetscErrorCode ierr; 1696ace3abfcSBarry Smith PetscBool match; 16974b9ad928SBarry Smith 16984b9ad928SBarry Smith PetscFunctionBegin; 16990700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 17004482741eSBarry Smith PetscValidIntPointer(n,2); 170192bb6962SLisandro Dalcin if (is) PetscValidPointer(is,3); 1702251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)pc,PCASM,&match);CHKERRQ(ierr); 17032a808120SBarry Smith if (!match) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"PC is not a PCASM"); 17044b9ad928SBarry Smith if (n) *n = osm->n_local_true; 17054b9ad928SBarry Smith if (is) *is = osm->is; 17062b691e39SMatthew Knepley if (is_local) *is_local = osm->is_local; 17074b9ad928SBarry Smith PetscFunctionReturn(0); 17084b9ad928SBarry Smith } 17094b9ad928SBarry Smith 17104b9ad928SBarry Smith /*@C 17114b9ad928SBarry Smith PCASMGetLocalSubmatrices - Gets the local submatrices (for this processor 17124b9ad928SBarry Smith only) for the additive Schwarz preconditioner. 17134b9ad928SBarry Smith 1714ad4df100SBarry Smith Not Collective 17154b9ad928SBarry Smith 17164b9ad928SBarry Smith Input Parameter: 17174b9ad928SBarry Smith . pc - the preconditioner context 17184b9ad928SBarry Smith 17194b9ad928SBarry Smith Output Parameters: 17204b9ad928SBarry Smith + n - the number of matrices for this processor (default value = 1) 17214b9ad928SBarry Smith - mat - the matrices 17224b9ad928SBarry Smith 17234b9ad928SBarry Smith Level: advanced 17244b9ad928SBarry Smith 172595452b02SPatrick Sanan Notes: 172634a84908Sprj- Call after PCSetUp() (or KSPSetUp()) but before PCApply() and before PCSetUpOnBlocks()) 1727cf739d55SBarry Smith 172834a84908Sprj- Usually one would use PCSetModifySubMatrices() to change the submatrices in building the preconditioner. 1729cf739d55SBarry Smith 17304b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(), 173134a84908Sprj- PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubdomains(), PCSetModifySubMatrices() 17324b9ad928SBarry Smith @*/ 17337087cfbeSBarry Smith PetscErrorCode PCASMGetLocalSubmatrices(PC pc,PetscInt *n,Mat *mat[]) 17344b9ad928SBarry Smith { 17354b9ad928SBarry Smith PC_ASM *osm; 173692bb6962SLisandro Dalcin PetscErrorCode ierr; 1737ace3abfcSBarry Smith PetscBool match; 17384b9ad928SBarry Smith 17394b9ad928SBarry Smith PetscFunctionBegin; 17400700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 174192bb6962SLisandro Dalcin PetscValidIntPointer(n,2); 174292bb6962SLisandro Dalcin if (mat) PetscValidPointer(mat,3); 174334a84908Sprj- if (!pc->setupcalled) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call after KSPSetUp() or PCSetUp()."); 1744251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)pc,PCASM,&match);CHKERRQ(ierr); 174592bb6962SLisandro Dalcin if (!match) { 174692bb6962SLisandro Dalcin if (n) *n = 0; 17470298fd71SBarry Smith if (mat) *mat = NULL; 174892bb6962SLisandro Dalcin } else { 17494b9ad928SBarry Smith osm = (PC_ASM*)pc->data; 17504b9ad928SBarry Smith if (n) *n = osm->n_local_true; 17514b9ad928SBarry Smith if (mat) *mat = osm->pmat; 175292bb6962SLisandro Dalcin } 17534b9ad928SBarry Smith PetscFunctionReturn(0); 17544b9ad928SBarry Smith } 1755d709ea83SDmitry Karpeev 1756d709ea83SDmitry Karpeev /*@ 1757d709ea83SDmitry Karpeev PCASMSetDMSubdomains - Indicates whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible. 1758f1ee410cSBarry Smith 1759d709ea83SDmitry Karpeev Logically Collective 1760d709ea83SDmitry Karpeev 1761d709ea83SDmitry Karpeev Input Parameter: 1762d709ea83SDmitry Karpeev + pc - the preconditioner 1763d709ea83SDmitry Karpeev - flg - boolean indicating whether to use subdomains defined by the DM 1764d709ea83SDmitry Karpeev 1765d709ea83SDmitry Karpeev Options Database Key: 1766d709ea83SDmitry Karpeev . -pc_asm_dm_subdomains 1767d709ea83SDmitry Karpeev 1768d709ea83SDmitry Karpeev Level: intermediate 1769d709ea83SDmitry Karpeev 1770d709ea83SDmitry Karpeev Notes: 1771d709ea83SDmitry Karpeev PCASMSetTotalSubdomains() and PCASMSetOverlap() take precedence over PCASMSetDMSubdomains(), 1772d709ea83SDmitry Karpeev so setting either of the first two effectively turns the latter off. 1773d709ea83SDmitry Karpeev 1774d709ea83SDmitry Karpeev .seealso: PCASMGetDMSubdomains(), PCASMSetTotalSubdomains(), PCASMSetOverlap() 1775d709ea83SDmitry Karpeev PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubdomains() 1776d709ea83SDmitry Karpeev @*/ 1777d709ea83SDmitry Karpeev PetscErrorCode PCASMSetDMSubdomains(PC pc,PetscBool flg) 1778d709ea83SDmitry Karpeev { 1779d709ea83SDmitry Karpeev PC_ASM *osm = (PC_ASM*)pc->data; 1780d709ea83SDmitry Karpeev PetscErrorCode ierr; 1781d709ea83SDmitry Karpeev PetscBool match; 1782d709ea83SDmitry Karpeev 1783d709ea83SDmitry Karpeev PetscFunctionBegin; 1784d709ea83SDmitry Karpeev PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1785d709ea83SDmitry Karpeev PetscValidLogicalCollectiveBool(pc,flg,2); 1786d709ea83SDmitry Karpeev if (pc->setupcalled) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Not for a setup PC."); 1787d709ea83SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)pc,PCASM,&match);CHKERRQ(ierr); 1788d709ea83SDmitry Karpeev if (match) { 1789d709ea83SDmitry Karpeev osm->dm_subdomains = flg; 1790d709ea83SDmitry Karpeev } 1791d709ea83SDmitry Karpeev PetscFunctionReturn(0); 1792d709ea83SDmitry Karpeev } 1793d709ea83SDmitry Karpeev 1794d709ea83SDmitry Karpeev /*@ 1795d709ea83SDmitry Karpeev PCASMGetDMSubdomains - Returns flag indicating whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible. 1796d709ea83SDmitry Karpeev Not Collective 1797d709ea83SDmitry Karpeev 1798d709ea83SDmitry Karpeev Input Parameter: 1799d709ea83SDmitry Karpeev . pc - the preconditioner 1800d709ea83SDmitry Karpeev 1801d709ea83SDmitry Karpeev Output Parameter: 1802d709ea83SDmitry Karpeev . flg - boolean indicating whether to use subdomains defined by the DM 1803d709ea83SDmitry Karpeev 1804d709ea83SDmitry Karpeev Level: intermediate 1805d709ea83SDmitry Karpeev 1806d709ea83SDmitry Karpeev .seealso: PCASMSetDMSubdomains(), PCASMSetTotalSubdomains(), PCASMSetOverlap() 1807d709ea83SDmitry Karpeev PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubdomains() 1808d709ea83SDmitry Karpeev @*/ 1809d709ea83SDmitry Karpeev PetscErrorCode PCASMGetDMSubdomains(PC pc,PetscBool* flg) 1810d709ea83SDmitry Karpeev { 1811d709ea83SDmitry Karpeev PC_ASM *osm = (PC_ASM*)pc->data; 1812d709ea83SDmitry Karpeev PetscErrorCode ierr; 1813d709ea83SDmitry Karpeev PetscBool match; 1814d709ea83SDmitry Karpeev 1815d709ea83SDmitry Karpeev PetscFunctionBegin; 1816d709ea83SDmitry Karpeev PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1817534a8f05SLisandro Dalcin PetscValidBoolPointer(flg,2); 1818d709ea83SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)pc,PCASM,&match);CHKERRQ(ierr); 1819d709ea83SDmitry Karpeev if (match) { 1820d709ea83SDmitry Karpeev if (flg) *flg = osm->dm_subdomains; 1821d709ea83SDmitry Karpeev } 1822d709ea83SDmitry Karpeev PetscFunctionReturn(0); 1823d709ea83SDmitry Karpeev } 182480ec0b7dSPatrick Sanan 182580ec0b7dSPatrick Sanan /*@ 182680ec0b7dSPatrick Sanan PCASMGetSubMatType - Gets the matrix type used for ASM subsolves, as a string. 182780ec0b7dSPatrick Sanan 182880ec0b7dSPatrick Sanan Not Collective 182980ec0b7dSPatrick Sanan 183080ec0b7dSPatrick Sanan Input Parameter: 183180ec0b7dSPatrick Sanan . pc - the PC 183280ec0b7dSPatrick Sanan 183380ec0b7dSPatrick Sanan Output Parameter: 1834f1ee410cSBarry Smith . -pc_asm_sub_mat_type - name of matrix type 183580ec0b7dSPatrick Sanan 183680ec0b7dSPatrick Sanan Level: advanced 183780ec0b7dSPatrick Sanan 183880ec0b7dSPatrick Sanan .seealso: PCASMSetSubMatType(), PCASM, PCSetType(), VecSetType(), MatType, Mat 183980ec0b7dSPatrick Sanan @*/ 184080ec0b7dSPatrick Sanan PetscErrorCode PCASMGetSubMatType(PC pc,MatType *sub_mat_type) { 184180ec0b7dSPatrick Sanan PetscErrorCode ierr; 184280ec0b7dSPatrick Sanan 184380ec0b7dSPatrick Sanan ierr = PetscTryMethod(pc,"PCASMGetSubMatType_C",(PC,MatType*),(pc,sub_mat_type));CHKERRQ(ierr); 184480ec0b7dSPatrick Sanan PetscFunctionReturn(0); 184580ec0b7dSPatrick Sanan } 184680ec0b7dSPatrick Sanan 184780ec0b7dSPatrick Sanan /*@ 184880ec0b7dSPatrick Sanan PCASMSetSubMatType - Set the type of matrix used for ASM subsolves 184980ec0b7dSPatrick Sanan 185080ec0b7dSPatrick Sanan Collective on Mat 185180ec0b7dSPatrick Sanan 185280ec0b7dSPatrick Sanan Input Parameters: 185380ec0b7dSPatrick Sanan + pc - the PC object 185480ec0b7dSPatrick Sanan - sub_mat_type - matrix type 185580ec0b7dSPatrick Sanan 185680ec0b7dSPatrick Sanan Options Database Key: 185780ec0b7dSPatrick 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. 185880ec0b7dSPatrick Sanan 185980ec0b7dSPatrick Sanan Notes: 186080ec0b7dSPatrick Sanan See "${PETSC_DIR}/include/petscmat.h" for available types 186180ec0b7dSPatrick Sanan 186280ec0b7dSPatrick Sanan Level: advanced 186380ec0b7dSPatrick Sanan 186480ec0b7dSPatrick Sanan .seealso: PCASMGetSubMatType(), PCASM, PCSetType(), VecSetType(), MatType, Mat 186580ec0b7dSPatrick Sanan @*/ 186680ec0b7dSPatrick Sanan PetscErrorCode PCASMSetSubMatType(PC pc,MatType sub_mat_type) 186780ec0b7dSPatrick Sanan { 186880ec0b7dSPatrick Sanan PetscErrorCode ierr; 186980ec0b7dSPatrick Sanan 187080ec0b7dSPatrick Sanan ierr = PetscTryMethod(pc,"PCASMSetSubMatType_C",(PC,MatType),(pc,sub_mat_type));CHKERRQ(ierr); 187180ec0b7dSPatrick Sanan PetscFunctionReturn(0); 187280ec0b7dSPatrick Sanan } 1873