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); 2821dd8081eSeaulisa 2834b9ad928SBarry Smith scall = MAT_INITIAL_MATRIX; 2844b9ad928SBarry Smith } else { 2854b9ad928SBarry Smith /* 2864b9ad928SBarry Smith Destroy the blocks from the previous iteration 2874b9ad928SBarry Smith */ 288e09e08ccSBarry Smith if (pc->flag == DIFFERENT_NONZERO_PATTERN) { 2894b9ad928SBarry Smith ierr = MatDestroyMatrices(osm->n_local_true,&osm->pmat);CHKERRQ(ierr); 2904b9ad928SBarry Smith scall = MAT_INITIAL_MATRIX; 2914b9ad928SBarry Smith } 2924b9ad928SBarry Smith } 2934b9ad928SBarry Smith 294*b58cb649SBarry Smith /* Destroy previous submatrices of a different type than pc->pmat since MAT_REUSE_MATRIX won't work in that case */ 295*b58cb649SBarry Smith if ((scall == MAT_REUSE_MATRIX) && osm->sub_mat_type) { 296*b58cb649SBarry Smith if (osm->n_local_true > 0) { 297*b58cb649SBarry Smith ierr = MatDestroySubMatrices(osm->n_local_true,&osm->pmat);CHKERRQ(ierr); 298*b58cb649SBarry Smith } 299*b58cb649SBarry Smith scall = MAT_INITIAL_MATRIX; 300*b58cb649SBarry Smith } 301*b58cb649SBarry Smith 30292bb6962SLisandro Dalcin /* 30392bb6962SLisandro Dalcin Extract out the submatrices 30492bb6962SLisandro Dalcin */ 3057dae84e0SHong Zhang ierr = MatCreateSubMatrices(pc->pmat,osm->n_local_true,osm->is,osm->is,scall,&osm->pmat);CHKERRQ(ierr); 30692bb6962SLisandro Dalcin if (scall == MAT_INITIAL_MATRIX) { 30792bb6962SLisandro Dalcin ierr = PetscObjectGetOptionsPrefix((PetscObject)pc->pmat,&pprefix);CHKERRQ(ierr); 30892bb6962SLisandro Dalcin for (i=0; i<osm->n_local_true; i++) { 3093bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)osm->pmat[i]);CHKERRQ(ierr); 31078904715SLisandro Dalcin ierr = PetscObjectSetOptionsPrefix((PetscObject)osm->pmat[i],pprefix);CHKERRQ(ierr); 31192bb6962SLisandro Dalcin } 31292bb6962SLisandro Dalcin } 31380ec0b7dSPatrick Sanan 31480ec0b7dSPatrick Sanan /* Convert the types of the submatrices (if needbe) */ 31580ec0b7dSPatrick Sanan if (osm->sub_mat_type) { 31680ec0b7dSPatrick Sanan for (i=0; i<osm->n_local_true; i++) { 31780ec0b7dSPatrick Sanan ierr = MatConvert(osm->pmat[i],osm->sub_mat_type,MAT_INPLACE_MATRIX,&(osm->pmat[i]));CHKERRQ(ierr); 31880ec0b7dSPatrick Sanan } 31980ec0b7dSPatrick Sanan } 32080ec0b7dSPatrick Sanan 32180ec0b7dSPatrick Sanan if (!pc->setupcalled) { 32256ea09ceSStefano Zampini VecType vtype; 32356ea09ceSStefano Zampini 32480ec0b7dSPatrick Sanan /* Create the local work vectors (from the local matrices) and scatter contexts */ 3250a545947SLisandro Dalcin ierr = MatCreateVecs(pc->pmat,&vec,NULL);CHKERRQ(ierr); 3265b3c0d42Seaulisa 3271dd8081eSeaulisa 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()"); 3281dd8081eSeaulisa if (osm->is_local && osm->type == PC_ASM_RESTRICT && osm->loctype == PC_COMPOSITE_ADDITIVE) { 3291dd8081eSeaulisa ierr = PetscMalloc1(osm->n_local_true,&osm->lprolongation);CHKERRQ(ierr); 3301dd8081eSeaulisa } 3311dd8081eSeaulisa ierr = PetscMalloc1(osm->n_local_true,&osm->lrestriction);CHKERRQ(ierr); 3321dd8081eSeaulisa ierr = PetscMalloc1(osm->n_local_true,&osm->x);CHKERRQ(ierr); 3331dd8081eSeaulisa ierr = PetscMalloc1(osm->n_local_true,&osm->y);CHKERRQ(ierr); 334347574c9Seaulisa 335347574c9Seaulisa ierr = ISGetLocalSize(osm->lis,&m);CHKERRQ(ierr); 336347574c9Seaulisa ierr = ISCreateStride(PETSC_COMM_SELF,m,0,1,&isl);CHKERRQ(ierr); 33756ea09ceSStefano Zampini ierr = MatGetVecType(osm->pmat[0],&vtype);CHKERRQ(ierr); 33856ea09ceSStefano Zampini ierr = VecCreate(PETSC_COMM_SELF,&osm->lx);CHKERRQ(ierr); 33956ea09ceSStefano Zampini ierr = VecSetSizes(osm->lx,m,m);CHKERRQ(ierr); 34056ea09ceSStefano Zampini ierr = VecSetType(osm->lx,vtype);CHKERRQ(ierr); 34156ea09ceSStefano Zampini ierr = VecDuplicate(osm->lx, &osm->ly);CHKERRQ(ierr); 3429448b7f1SJunchao Zhang ierr = VecScatterCreate(vec,osm->lis,osm->lx,isl,&osm->restriction);CHKERRQ(ierr); 343347574c9Seaulisa ierr = ISDestroy(&isl);CHKERRQ(ierr); 344347574c9Seaulisa 34580ec0b7dSPatrick Sanan for (i=0; i<osm->n_local_true; ++i) { 3465b3c0d42Seaulisa ISLocalToGlobalMapping ltog; 3475b3c0d42Seaulisa IS isll; 3485b3c0d42Seaulisa const PetscInt *idx_is; 3495b3c0d42Seaulisa PetscInt *idx_lis,nout; 3505b3c0d42Seaulisa 35155817e79SBarry Smith ierr = ISGetLocalSize(osm->is[i],&m);CHKERRQ(ierr); 35255817e79SBarry Smith ierr = MatCreateVecs(osm->pmat[i],&osm->x[i],NULL);CHKERRQ(ierr); 35355817e79SBarry Smith ierr = VecDuplicate(osm->x[i],&osm->y[i]);CHKERRQ(ierr); 35455817e79SBarry Smith 355b0de9f23SBarry Smith /* generate a scatter from ly to y[i] picking all the overlapping is[i] entries */ 3565b3c0d42Seaulisa ierr = ISLocalToGlobalMappingCreateIS(osm->lis,<og);CHKERRQ(ierr); 3575b3c0d42Seaulisa ierr = ISGetLocalSize(osm->is[i],&m);CHKERRQ(ierr); 3585b3c0d42Seaulisa ierr = ISGetIndices(osm->is[i], &idx_is);CHKERRQ(ierr); 3595b3c0d42Seaulisa ierr = PetscMalloc1(m,&idx_lis);CHKERRQ(ierr); 3605b3c0d42Seaulisa ierr = ISGlobalToLocalMappingApply(ltog,IS_GTOLM_DROP,m,idx_is,&nout,idx_lis);CHKERRQ(ierr); 3615b3c0d42Seaulisa if (nout != m) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"is not a subset of lis"); 3625b3c0d42Seaulisa ierr = ISRestoreIndices(osm->is[i], &idx_is);CHKERRQ(ierr); 3635b3c0d42Seaulisa ierr = ISCreateGeneral(PETSC_COMM_SELF,m,idx_lis,PETSC_OWN_POINTER,&isll);CHKERRQ(ierr); 364d8b3b5e3Seaulisa ierr = ISLocalToGlobalMappingDestroy(<og);CHKERRQ(ierr); 3655b3c0d42Seaulisa ierr = ISCreateStride(PETSC_COMM_SELF,m,0,1,&isl);CHKERRQ(ierr); 3669448b7f1SJunchao Zhang ierr = VecScatterCreate(osm->ly,isll,osm->y[i],isl,&osm->lrestriction[i]);CHKERRQ(ierr); 3675b3c0d42Seaulisa ierr = ISDestroy(&isll);CHKERRQ(ierr); 3685b3c0d42Seaulisa ierr = ISDestroy(&isl);CHKERRQ(ierr); 369910cf402Sprj- if (osm->lprolongation) { /* generate a scatter from y[i] to ly picking only the the non-overlapping is_local[i] entries */ 370d8b3b5e3Seaulisa ISLocalToGlobalMapping ltog; 371d8b3b5e3Seaulisa IS isll,isll_local; 372d8b3b5e3Seaulisa const PetscInt *idx_local; 373d8b3b5e3Seaulisa PetscInt *idx1, *idx2, nout; 374d8b3b5e3Seaulisa 375d8b3b5e3Seaulisa ierr = ISGetLocalSize(osm->is_local[i],&m_local);CHKERRQ(ierr); 376d8b3b5e3Seaulisa ierr = ISGetIndices(osm->is_local[i], &idx_local);CHKERRQ(ierr); 377d8b3b5e3Seaulisa 378d8b3b5e3Seaulisa ierr = ISLocalToGlobalMappingCreateIS(osm->is[i],<og);CHKERRQ(ierr); 379d8b3b5e3Seaulisa ierr = PetscMalloc1(m_local,&idx1);CHKERRQ(ierr); 380d8b3b5e3Seaulisa ierr = ISGlobalToLocalMappingApply(ltog,IS_GTOLM_DROP,m_local,idx_local,&nout,idx1);CHKERRQ(ierr); 381d8b3b5e3Seaulisa ierr = ISLocalToGlobalMappingDestroy(<og);CHKERRQ(ierr); 382d8b3b5e3Seaulisa if (nout != m_local) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"is_local not a subset of is"); 383d8b3b5e3Seaulisa ierr = ISCreateGeneral(PETSC_COMM_SELF,m_local,idx1,PETSC_OWN_POINTER,&isll);CHKERRQ(ierr); 384d8b3b5e3Seaulisa 385d8b3b5e3Seaulisa ierr = ISLocalToGlobalMappingCreateIS(osm->lis,<og);CHKERRQ(ierr); 386d8b3b5e3Seaulisa ierr = PetscMalloc1(m_local,&idx2);CHKERRQ(ierr); 387d8b3b5e3Seaulisa ierr = ISGlobalToLocalMappingApply(ltog,IS_GTOLM_DROP,m_local,idx_local,&nout,idx2);CHKERRQ(ierr); 388d8b3b5e3Seaulisa ierr = ISLocalToGlobalMappingDestroy(<og);CHKERRQ(ierr); 389d8b3b5e3Seaulisa if (nout != m_local) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_PLIB,"is_local not a subset of lis"); 390d8b3b5e3Seaulisa ierr = ISCreateGeneral(PETSC_COMM_SELF,m_local,idx2,PETSC_OWN_POINTER,&isll_local);CHKERRQ(ierr); 391d8b3b5e3Seaulisa 392d8b3b5e3Seaulisa ierr = ISRestoreIndices(osm->is_local[i], &idx_local);CHKERRQ(ierr); 3939448b7f1SJunchao Zhang ierr = VecScatterCreate(osm->y[i],isll,osm->ly,isll_local,&osm->lprolongation[i]);CHKERRQ(ierr); 394d8b3b5e3Seaulisa 395d8b3b5e3Seaulisa ierr = ISDestroy(&isll);CHKERRQ(ierr); 396d8b3b5e3Seaulisa ierr = ISDestroy(&isll_local);CHKERRQ(ierr); 397d8b3b5e3Seaulisa } 39880ec0b7dSPatrick Sanan } 39980ec0b7dSPatrick Sanan ierr = VecDestroy(&vec);CHKERRQ(ierr); 40080ec0b7dSPatrick Sanan } 40180ec0b7dSPatrick Sanan 402fb745f2cSMatthew G. Knepley if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE) { 403235cc4ceSMatthew G. Knepley IS *cis; 404235cc4ceSMatthew G. Knepley PetscInt c; 405235cc4ceSMatthew G. Knepley 406235cc4ceSMatthew G. Knepley ierr = PetscMalloc1(osm->n_local_true, &cis);CHKERRQ(ierr); 407235cc4ceSMatthew G. Knepley for (c = 0; c < osm->n_local_true; ++c) cis[c] = osm->lis; 4087dae84e0SHong Zhang ierr = MatCreateSubMatrices(pc->pmat, osm->n_local_true, osm->is, cis, scall, &osm->lmats);CHKERRQ(ierr); 409235cc4ceSMatthew G. Knepley ierr = PetscFree(cis);CHKERRQ(ierr); 410fb745f2cSMatthew G. Knepley } 4114b9ad928SBarry Smith 4124b9ad928SBarry Smith /* Return control to the user so that the submatrices can be modified (e.g., to apply 4134b9ad928SBarry Smith different boundary conditions for the submatrices than for the global problem) */ 41492bb6962SLisandro Dalcin ierr = PCModifySubMatrices(pc,osm->n_local_true,osm->is,osm->is,osm->pmat,pc->modifysubmatricesP);CHKERRQ(ierr); 4154b9ad928SBarry Smith 41692bb6962SLisandro Dalcin /* 41792bb6962SLisandro Dalcin Loop over subdomains putting them into local ksp 41892bb6962SLisandro Dalcin */ 41992bb6962SLisandro Dalcin for (i=0; i<osm->n_local_true; i++) { 42023ee1639SBarry Smith ierr = KSPSetOperators(osm->ksp[i],osm->pmat[i],osm->pmat[i]);CHKERRQ(ierr); 42192bb6962SLisandro Dalcin if (!pc->setupcalled) { 422bf108f30SBarry Smith ierr = KSPSetFromOptions(osm->ksp[i]);CHKERRQ(ierr); 4234b9ad928SBarry Smith } 42492bb6962SLisandro Dalcin } 4254b9ad928SBarry Smith PetscFunctionReturn(0); 4264b9ad928SBarry Smith } 4274b9ad928SBarry Smith 4286849ba73SBarry Smith static PetscErrorCode PCSetUpOnBlocks_ASM(PC pc) 4294b9ad928SBarry Smith { 4304b9ad928SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 4316849ba73SBarry Smith PetscErrorCode ierr; 43213f74950SBarry Smith PetscInt i; 433539c167fSBarry Smith KSPConvergedReason reason; 4344b9ad928SBarry Smith 4354b9ad928SBarry Smith PetscFunctionBegin; 4364b9ad928SBarry Smith for (i=0; i<osm->n_local_true; i++) { 4374b9ad928SBarry Smith ierr = KSPSetUp(osm->ksp[i]);CHKERRQ(ierr); 438539c167fSBarry Smith ierr = KSPGetConvergedReason(osm->ksp[i],&reason);CHKERRQ(ierr); 439c0decd05SBarry Smith if (reason == KSP_DIVERGED_PC_FAILED) { 440261222b2SHong Zhang pc->failedreason = PC_SUBPC_ERROR; 441e0eafd54SHong Zhang } 4424b9ad928SBarry Smith } 4434b9ad928SBarry Smith PetscFunctionReturn(0); 4444b9ad928SBarry Smith } 4454b9ad928SBarry Smith 4466849ba73SBarry Smith static PetscErrorCode PCApply_ASM(PC pc,Vec x,Vec y) 4474b9ad928SBarry Smith { 4484b9ad928SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 4496849ba73SBarry Smith PetscErrorCode ierr; 4501dd8081eSeaulisa PetscInt i,n_local_true = osm->n_local_true; 4514b9ad928SBarry Smith ScatterMode forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE; 4524b9ad928SBarry Smith 4534b9ad928SBarry Smith PetscFunctionBegin; 4544b9ad928SBarry Smith /* 45548e38a8aSPierre Jolivet support for limiting the restriction or interpolation to only local 4564b9ad928SBarry Smith subdomain values (leaving the other values 0). 4574b9ad928SBarry Smith */ 4584b9ad928SBarry Smith if (!(osm->type & PC_ASM_RESTRICT)) { 4594b9ad928SBarry Smith forward = SCATTER_FORWARD_LOCAL; 4604b9ad928SBarry Smith /* have to zero the work RHS since scatter may leave some slots empty */ 4611dd8081eSeaulisa ierr = VecSet(osm->lx, 0.0);CHKERRQ(ierr); 4624b9ad928SBarry Smith } 463347574c9Seaulisa if (!(osm->type & PC_ASM_INTERPOLATE)) { 464347574c9Seaulisa reverse = SCATTER_REVERSE_LOCAL; 465347574c9Seaulisa } 4664b9ad928SBarry Smith 467347574c9Seaulisa if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE || osm->loctype == PC_COMPOSITE_ADDITIVE) { 468b0de9f23SBarry Smith /* zero the global and the local solutions */ 469910cf402Sprj- ierr = VecSet(y, 0.0);CHKERRQ(ierr); 4705b3c0d42Seaulisa ierr = VecSet(osm->ly, 0.0);CHKERRQ(ierr); 471347574c9Seaulisa 47248e38a8aSPierre Jolivet /* copy the global RHS to local RHS including the ghost nodes */ 4731dd8081eSeaulisa ierr = VecScatterBegin(osm->restriction, x, osm->lx, INSERT_VALUES, forward);CHKERRQ(ierr); 4741dd8081eSeaulisa ierr = VecScatterEnd(osm->restriction, x, osm->lx, INSERT_VALUES, forward);CHKERRQ(ierr); 475347574c9Seaulisa 47648e38a8aSPierre Jolivet /* restrict local RHS to the overlapping 0-block RHS */ 4771dd8081eSeaulisa ierr = VecScatterBegin(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward);CHKERRQ(ierr); 4781dd8081eSeaulisa ierr = VecScatterEnd(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward);CHKERRQ(ierr); 479d8b3b5e3Seaulisa 48012cd4985SMatthew G. Knepley /* do the local solves */ 48112cd4985SMatthew G. Knepley for (i = 0; i < n_local_true; ++i) { 482347574c9Seaulisa 483b0de9f23SBarry Smith /* solve the overlapping i-block */ 484fa2bb9feSLisandro Dalcin ierr = PetscLogEventBegin(PC_ApplyOnBlocks, osm->ksp[i], osm->x[i], osm->y[i],0);CHKERRQ(ierr); 48512cd4985SMatthew G. Knepley ierr = KSPSolve(osm->ksp[i], osm->x[i], osm->y[i]);CHKERRQ(ierr); 486c0decd05SBarry Smith ierr = KSPCheckSolve(osm->ksp[i], pc, osm->y[i]);CHKERRQ(ierr); 487fa2bb9feSLisandro Dalcin ierr = PetscLogEventEnd(PC_ApplyOnBlocks, osm->ksp[i], osm->x[i], osm->y[i], 0);CHKERRQ(ierr); 488d8b3b5e3Seaulisa 489910cf402Sprj- if (osm->lprolongation) { /* interpolate the non-overlapping i-block solution to the local solution (only for restrictive additive) */ 4901dd8081eSeaulisa ierr = VecScatterBegin(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward);CHKERRQ(ierr); 4911dd8081eSeaulisa ierr = VecScatterEnd(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward);CHKERRQ(ierr); 49248e38a8aSPierre Jolivet } else { /* interpolate the overlapping i-block solution to the local solution */ 4931dd8081eSeaulisa ierr = VecScatterBegin(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse);CHKERRQ(ierr); 4941dd8081eSeaulisa ierr = VecScatterEnd(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse);CHKERRQ(ierr); 495d8b3b5e3Seaulisa } 496347574c9Seaulisa 497347574c9Seaulisa if (i < n_local_true-1) { 49848e38a8aSPierre Jolivet /* restrict local RHS to the overlapping (i+1)-block RHS */ 4991dd8081eSeaulisa ierr = VecScatterBegin(osm->lrestriction[i+1], osm->lx, osm->x[i+1], INSERT_VALUES, forward);CHKERRQ(ierr); 5001dd8081eSeaulisa ierr = VecScatterEnd(osm->lrestriction[i+1], osm->lx, osm->x[i+1], INSERT_VALUES, forward);CHKERRQ(ierr); 501347574c9Seaulisa 502347574c9Seaulisa if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE) { 5038966356dSPierre Jolivet /* update the overlapping (i+1)-block RHS using the current local solution */ 504347574c9Seaulisa ierr = MatMult(osm->lmats[i+1], osm->ly, osm->y[i+1]);CHKERRQ(ierr); 505347574c9Seaulisa ierr = VecAXPBY(osm->x[i+1],-1.,1., osm->y[i+1]);CHKERRQ(ierr); 5067c3d802fSMatthew G. Knepley } 50712cd4985SMatthew G. Knepley } 50812cd4985SMatthew G. Knepley } 50948e38a8aSPierre Jolivet /* add the local solution to the global solution including the ghost nodes */ 5101dd8081eSeaulisa ierr = VecScatterBegin(osm->restriction, osm->ly, y, ADD_VALUES, reverse);CHKERRQ(ierr); 5111dd8081eSeaulisa ierr = VecScatterEnd(osm->restriction, osm->ly, y, ADD_VALUES, reverse);CHKERRQ(ierr); 51289c010cfSBarry Smith } else SETERRQ1(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Invalid local composition type: %s", PCCompositeTypes[osm->loctype]); 5134b9ad928SBarry Smith PetscFunctionReturn(0); 5144b9ad928SBarry Smith } 5154b9ad928SBarry Smith 51648e38a8aSPierre Jolivet static PetscErrorCode PCMatApply_ASM(PC pc,Mat X,Mat Y) 51748e38a8aSPierre Jolivet { 51848e38a8aSPierre Jolivet PC_ASM *osm = (PC_ASM*)pc->data; 51948e38a8aSPierre Jolivet Mat Z,W; 52048e38a8aSPierre Jolivet Vec x; 52148e38a8aSPierre Jolivet PetscInt i,m,N; 52248e38a8aSPierre Jolivet ScatterMode forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE; 52348e38a8aSPierre Jolivet PetscErrorCode ierr; 52448e38a8aSPierre Jolivet 52548e38a8aSPierre Jolivet PetscFunctionBegin; 52648e38a8aSPierre Jolivet if (osm->n_local_true > 1) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Not yet implemented"); 52748e38a8aSPierre Jolivet /* 52848e38a8aSPierre Jolivet support for limiting the restriction or interpolation to only local 52948e38a8aSPierre Jolivet subdomain values (leaving the other values 0). 53048e38a8aSPierre Jolivet */ 53148e38a8aSPierre Jolivet if (!(osm->type & PC_ASM_RESTRICT)) { 53248e38a8aSPierre Jolivet forward = SCATTER_FORWARD_LOCAL; 53348e38a8aSPierre Jolivet /* have to zero the work RHS since scatter may leave some slots empty */ 53448e38a8aSPierre Jolivet ierr = VecSet(osm->lx, 0.0);CHKERRQ(ierr); 53548e38a8aSPierre Jolivet } 53648e38a8aSPierre Jolivet if (!(osm->type & PC_ASM_INTERPOLATE)) { 53748e38a8aSPierre Jolivet reverse = SCATTER_REVERSE_LOCAL; 53848e38a8aSPierre Jolivet } 53948e38a8aSPierre Jolivet ierr = VecGetLocalSize(osm->x[0], &m);CHKERRQ(ierr); 54048e38a8aSPierre Jolivet ierr = MatGetSize(X, NULL, &N);CHKERRQ(ierr); 54148e38a8aSPierre Jolivet ierr = MatCreateSeqDense(PETSC_COMM_SELF, m, N, NULL, &Z);CHKERRQ(ierr); 54248e38a8aSPierre Jolivet if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE || osm->loctype == PC_COMPOSITE_ADDITIVE) { 54348e38a8aSPierre Jolivet /* zero the global and the local solutions */ 54448e38a8aSPierre Jolivet ierr = MatZeroEntries(Y);CHKERRQ(ierr); 54548e38a8aSPierre Jolivet ierr = VecSet(osm->ly, 0.0);CHKERRQ(ierr); 54648e38a8aSPierre Jolivet 54748e38a8aSPierre Jolivet for (i = 0; i < N; ++i) { 54848e38a8aSPierre Jolivet ierr = MatDenseGetColumnVecRead(X, i, &x); 54948e38a8aSPierre Jolivet /* copy the global RHS to local RHS including the ghost nodes */ 55048e38a8aSPierre Jolivet ierr = VecScatterBegin(osm->restriction, x, osm->lx, INSERT_VALUES, forward);CHKERRQ(ierr); 55148e38a8aSPierre Jolivet ierr = VecScatterEnd(osm->restriction, x, osm->lx, INSERT_VALUES, forward);CHKERRQ(ierr); 55248e38a8aSPierre Jolivet ierr = MatDenseRestoreColumnVecRead(X, i, &x); 55348e38a8aSPierre Jolivet 55448e38a8aSPierre Jolivet ierr = MatDenseGetColumnVecWrite(Z, i, &x); 55548e38a8aSPierre Jolivet /* restrict local RHS to the overlapping 0-block RHS */ 55648e38a8aSPierre Jolivet ierr = VecScatterBegin(osm->lrestriction[0], osm->lx, x, INSERT_VALUES, forward);CHKERRQ(ierr); 55748e38a8aSPierre Jolivet ierr = VecScatterEnd(osm->lrestriction[0], osm->lx, x, INSERT_VALUES, forward);CHKERRQ(ierr); 55848e38a8aSPierre Jolivet ierr = MatDenseRestoreColumnVecWrite(Z, i, &x); 55948e38a8aSPierre Jolivet } 56048e38a8aSPierre Jolivet ierr = MatCreateSeqDense(PETSC_COMM_SELF, m, N, NULL, &W);CHKERRQ(ierr); 56148e38a8aSPierre Jolivet /* solve the overlapping 0-block */ 56248e38a8aSPierre Jolivet ierr = PetscLogEventBegin(PC_ApplyOnBlocks, osm->ksp[0], Z, W, 0);CHKERRQ(ierr); 56348e38a8aSPierre Jolivet ierr = KSPMatSolve(osm->ksp[0], Z, W);CHKERRQ(ierr); 56448e38a8aSPierre Jolivet ierr = KSPCheckSolve(osm->ksp[0], pc, NULL);CHKERRQ(ierr); 56548e38a8aSPierre Jolivet ierr = PetscLogEventEnd(PC_ApplyOnBlocks, osm->ksp[0], Z, W,0);CHKERRQ(ierr); 56648e38a8aSPierre Jolivet ierr = MatDestroy(&Z);CHKERRQ(ierr); 56748e38a8aSPierre Jolivet 56848e38a8aSPierre Jolivet for (i = 0; i < N; ++i) { 56948e38a8aSPierre Jolivet ierr = VecSet(osm->ly, 0.0);CHKERRQ(ierr); 57048e38a8aSPierre Jolivet ierr = MatDenseGetColumnVecRead(W, i, &x); 57148e38a8aSPierre Jolivet if (osm->lprolongation) { /* interpolate the non-overlapping 0-block solution to the local solution (only for restrictive additive) */ 57248e38a8aSPierre Jolivet ierr = VecScatterBegin(osm->lprolongation[0], x, osm->ly, ADD_VALUES, forward);CHKERRQ(ierr); 57348e38a8aSPierre Jolivet ierr = VecScatterEnd(osm->lprolongation[0], x, osm->ly, ADD_VALUES, forward);CHKERRQ(ierr); 57448e38a8aSPierre Jolivet } else { /* interpolate the overlapping 0-block solution to the local solution */ 57548e38a8aSPierre Jolivet ierr = VecScatterBegin(osm->lrestriction[0], x, osm->ly, ADD_VALUES, reverse);CHKERRQ(ierr); 57648e38a8aSPierre Jolivet ierr = VecScatterEnd(osm->lrestriction[0], x, osm->ly, ADD_VALUES, reverse);CHKERRQ(ierr); 57748e38a8aSPierre Jolivet } 57848e38a8aSPierre Jolivet ierr = MatDenseRestoreColumnVecRead(W, i, &x); 57948e38a8aSPierre Jolivet 58048e38a8aSPierre Jolivet ierr = MatDenseGetColumnVecWrite(Y, i, &x); 58148e38a8aSPierre Jolivet /* add the local solution to the global solution including the ghost nodes */ 58248e38a8aSPierre Jolivet ierr = VecScatterBegin(osm->restriction, osm->ly, x, ADD_VALUES, reverse);CHKERRQ(ierr); 58348e38a8aSPierre Jolivet ierr = VecScatterEnd(osm->restriction, osm->ly, x, ADD_VALUES, reverse);CHKERRQ(ierr); 58448e38a8aSPierre Jolivet ierr = MatDenseRestoreColumnVecWrite(Y, i, &x); 58548e38a8aSPierre Jolivet } 58648e38a8aSPierre Jolivet ierr = MatDestroy(&W);CHKERRQ(ierr); 58789c010cfSBarry Smith } else SETERRQ1(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Invalid local composition type: %s", PCCompositeTypes[osm->loctype]); 58848e38a8aSPierre Jolivet PetscFunctionReturn(0); 58948e38a8aSPierre Jolivet } 59048e38a8aSPierre Jolivet 5916849ba73SBarry Smith static PetscErrorCode PCApplyTranspose_ASM(PC pc,Vec x,Vec y) 5924b9ad928SBarry Smith { 5934b9ad928SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 5946849ba73SBarry Smith PetscErrorCode ierr; 5951dd8081eSeaulisa PetscInt i,n_local_true = osm->n_local_true; 5964b9ad928SBarry Smith ScatterMode forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE; 5974b9ad928SBarry Smith 5984b9ad928SBarry Smith PetscFunctionBegin; 5994b9ad928SBarry Smith /* 6004b9ad928SBarry Smith Support for limiting the restriction or interpolation to only local 6014b9ad928SBarry Smith subdomain values (leaving the other values 0). 6024b9ad928SBarry Smith 6034b9ad928SBarry Smith Note: these are reversed from the PCApply_ASM() because we are applying the 6044b9ad928SBarry Smith transpose of the three terms 6054b9ad928SBarry Smith */ 606d8b3b5e3Seaulisa 6074b9ad928SBarry Smith if (!(osm->type & PC_ASM_INTERPOLATE)) { 6084b9ad928SBarry Smith forward = SCATTER_FORWARD_LOCAL; 6094b9ad928SBarry Smith /* have to zero the work RHS since scatter may leave some slots empty */ 6101dd8081eSeaulisa ierr = VecSet(osm->lx, 0.0);CHKERRQ(ierr); 6114b9ad928SBarry Smith } 6122fa5cd67SKarl Rupp if (!(osm->type & PC_ASM_RESTRICT)) reverse = SCATTER_REVERSE_LOCAL; 6134b9ad928SBarry Smith 614b0de9f23SBarry Smith /* zero the global and the local solutions */ 615910cf402Sprj- ierr = VecSet(y, 0.0);CHKERRQ(ierr); 616d8b3b5e3Seaulisa ierr = VecSet(osm->ly, 0.0);CHKERRQ(ierr); 617d8b3b5e3Seaulisa 618b0de9f23SBarry Smith /* Copy the global RHS to local RHS including the ghost nodes */ 6191dd8081eSeaulisa ierr = VecScatterBegin(osm->restriction, x, osm->lx, INSERT_VALUES, forward);CHKERRQ(ierr); 6201dd8081eSeaulisa ierr = VecScatterEnd(osm->restriction, x, osm->lx, INSERT_VALUES, forward);CHKERRQ(ierr); 621d8b3b5e3Seaulisa 622b0de9f23SBarry Smith /* Restrict local RHS to the overlapping 0-block RHS */ 6231dd8081eSeaulisa ierr = VecScatterBegin(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward);CHKERRQ(ierr); 6241dd8081eSeaulisa ierr = VecScatterEnd(osm->lrestriction[0], osm->lx, osm->x[0], INSERT_VALUES, forward);CHKERRQ(ierr); 625d8b3b5e3Seaulisa 6264b9ad928SBarry Smith /* do the local solves */ 627d8b3b5e3Seaulisa for (i = 0; i < n_local_true; ++i) { 628d8b3b5e3Seaulisa 629b0de9f23SBarry Smith /* solve the overlapping i-block */ 630fa2bb9feSLisandro Dalcin ierr = PetscLogEventBegin(PC_ApplyOnBlocks,osm->ksp[i],osm->x[i],osm->y[i],0);CHKERRQ(ierr); 631e1bcd54cSBarry Smith ierr = KSPSolveTranspose(osm->ksp[i], osm->x[i], osm->y[i]);CHKERRQ(ierr); 632c0decd05SBarry Smith ierr = KSPCheckSolve(osm->ksp[i],pc,osm->y[i]);CHKERRQ(ierr); 633fa2bb9feSLisandro Dalcin ierr = PetscLogEventEnd(PC_ApplyOnBlocks,osm->ksp[i],osm->x[i],osm->y[i],0);CHKERRQ(ierr); 634d8b3b5e3Seaulisa 635910cf402Sprj- if (osm->lprolongation) { /* interpolate the non-overlapping i-block solution to the local solution */ 6361dd8081eSeaulisa ierr = VecScatterBegin(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward);CHKERRQ(ierr); 6371dd8081eSeaulisa ierr = VecScatterEnd(osm->lprolongation[i], osm->y[i], osm->ly, ADD_VALUES, forward);CHKERRQ(ierr); 638910cf402Sprj- } else { /* interpolate the overlapping i-block solution to the local solution */ 6391dd8081eSeaulisa ierr = VecScatterBegin(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse);CHKERRQ(ierr); 6401dd8081eSeaulisa ierr = VecScatterEnd(osm->lrestriction[i], osm->y[i], osm->ly, ADD_VALUES, reverse);CHKERRQ(ierr); 6414b9ad928SBarry Smith } 642d8b3b5e3Seaulisa 643d8b3b5e3Seaulisa if (i < n_local_true-1) { 644b0de9f23SBarry Smith /* Restrict local RHS to the overlapping (i+1)-block RHS */ 6451dd8081eSeaulisa ierr = VecScatterBegin(osm->lrestriction[i+1], osm->lx, osm->x[i+1], INSERT_VALUES, forward);CHKERRQ(ierr); 6461dd8081eSeaulisa ierr = VecScatterEnd(osm->lrestriction[i+1], osm->lx, osm->x[i+1], INSERT_VALUES, forward);CHKERRQ(ierr); 6474b9ad928SBarry Smith } 6484b9ad928SBarry Smith } 649b0de9f23SBarry Smith /* Add the local solution to the global solution including the ghost nodes */ 6501dd8081eSeaulisa ierr = VecScatterBegin(osm->restriction, osm->ly, y, ADD_VALUES, reverse);CHKERRQ(ierr); 6511dd8081eSeaulisa ierr = VecScatterEnd(osm->restriction, osm->ly, y, ADD_VALUES, reverse);CHKERRQ(ierr); 652d8b3b5e3Seaulisa 6534b9ad928SBarry Smith PetscFunctionReturn(0); 654d8b3b5e3Seaulisa 6554b9ad928SBarry Smith } 6564b9ad928SBarry Smith 657e91c6855SBarry Smith static PetscErrorCode PCReset_ASM(PC pc) 6584b9ad928SBarry Smith { 6594b9ad928SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 6606849ba73SBarry Smith PetscErrorCode ierr; 66113f74950SBarry Smith PetscInt i; 6624b9ad928SBarry Smith 6634b9ad928SBarry Smith PetscFunctionBegin; 66492bb6962SLisandro Dalcin if (osm->ksp) { 66592bb6962SLisandro Dalcin for (i=0; i<osm->n_local_true; i++) { 666e91c6855SBarry Smith ierr = KSPReset(osm->ksp[i]);CHKERRQ(ierr); 66792bb6962SLisandro Dalcin } 66892bb6962SLisandro Dalcin } 669e09e08ccSBarry Smith if (osm->pmat) { 67092bb6962SLisandro Dalcin if (osm->n_local_true > 0) { 67130a70a9aSHong Zhang ierr = MatDestroySubMatrices(osm->n_local_true,&osm->pmat);CHKERRQ(ierr); 67292bb6962SLisandro Dalcin } 67392bb6962SLisandro Dalcin } 6741dd8081eSeaulisa if (osm->lrestriction) { 6751dd8081eSeaulisa ierr = VecScatterDestroy(&osm->restriction);CHKERRQ(ierr); 6761dd8081eSeaulisa for (i=0; i<osm->n_local_true; i++) { 6771dd8081eSeaulisa ierr = VecScatterDestroy(&osm->lrestriction[i]);CHKERRQ(ierr); 6781dd8081eSeaulisa if (osm->lprolongation) {ierr = VecScatterDestroy(&osm->lprolongation[i]);CHKERRQ(ierr);} 679fcfd50ebSBarry Smith ierr = VecDestroy(&osm->x[i]);CHKERRQ(ierr); 680fcfd50ebSBarry Smith ierr = VecDestroy(&osm->y[i]);CHKERRQ(ierr); 6814b9ad928SBarry Smith } 6821dd8081eSeaulisa ierr = PetscFree(osm->lrestriction);CHKERRQ(ierr); 6831dd8081eSeaulisa if (osm->lprolongation) {ierr = PetscFree(osm->lprolongation);CHKERRQ(ierr);} 68405b42c5fSBarry Smith ierr = PetscFree(osm->x);CHKERRQ(ierr); 68578904715SLisandro Dalcin ierr = PetscFree(osm->y);CHKERRQ(ierr); 6861dd8081eSeaulisa 68792bb6962SLisandro Dalcin } 6882b691e39SMatthew Knepley ierr = PCASMDestroySubdomains(osm->n_local_true,osm->is,osm->is_local);CHKERRQ(ierr); 689fb745f2cSMatthew G. Knepley ierr = ISDestroy(&osm->lis);CHKERRQ(ierr); 690fb745f2cSMatthew G. Knepley ierr = VecDestroy(&osm->lx);CHKERRQ(ierr); 691fb745f2cSMatthew G. Knepley ierr = VecDestroy(&osm->ly);CHKERRQ(ierr); 692347574c9Seaulisa if (osm->loctype == PC_COMPOSITE_MULTIPLICATIVE) { 693347574c9Seaulisa ierr = MatDestroyMatrices(osm->n_local_true, &osm->lmats);CHKERRQ(ierr); 694fb745f2cSMatthew G. Knepley } 6952fa5cd67SKarl Rupp 69680ec0b7dSPatrick Sanan ierr = PetscFree(osm->sub_mat_type);CHKERRQ(ierr); 69780ec0b7dSPatrick Sanan 6980a545947SLisandro Dalcin osm->is = NULL; 6990a545947SLisandro Dalcin osm->is_local = NULL; 700e91c6855SBarry Smith PetscFunctionReturn(0); 701e91c6855SBarry Smith } 702e91c6855SBarry Smith 703e91c6855SBarry Smith static PetscErrorCode PCDestroy_ASM(PC pc) 704e91c6855SBarry Smith { 705e91c6855SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 706e91c6855SBarry Smith PetscErrorCode ierr; 707e91c6855SBarry Smith PetscInt i; 708e91c6855SBarry Smith 709e91c6855SBarry Smith PetscFunctionBegin; 710e91c6855SBarry Smith ierr = PCReset_ASM(pc);CHKERRQ(ierr); 711e91c6855SBarry Smith if (osm->ksp) { 712e91c6855SBarry Smith for (i=0; i<osm->n_local_true; i++) { 713fcfd50ebSBarry Smith ierr = KSPDestroy(&osm->ksp[i]);CHKERRQ(ierr); 714e91c6855SBarry Smith } 715e91c6855SBarry Smith ierr = PetscFree(osm->ksp);CHKERRQ(ierr); 716e91c6855SBarry Smith } 717e91c6855SBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 7184b9ad928SBarry Smith PetscFunctionReturn(0); 7194b9ad928SBarry Smith } 7204b9ad928SBarry Smith 7214416b707SBarry Smith static PetscErrorCode PCSetFromOptions_ASM(PetscOptionItems *PetscOptionsObject,PC pc) 7224b9ad928SBarry Smith { 7234b9ad928SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 7246849ba73SBarry Smith PetscErrorCode ierr; 7259dcbbd2bSBarry Smith PetscInt blocks,ovl; 72657501b6eSBarry Smith PetscBool flg; 72792bb6962SLisandro Dalcin PCASMType asmtype; 72812cd4985SMatthew G. Knepley PCCompositeType loctype; 72980ec0b7dSPatrick Sanan char sub_mat_type[256]; 7304b9ad928SBarry Smith 7314b9ad928SBarry Smith PetscFunctionBegin; 732e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"Additive Schwarz options");CHKERRQ(ierr); 733d709ea83SDmitry Karpeev ierr = PetscOptionsBool("-pc_asm_dm_subdomains","Use DMCreateDomainDecomposition() to define subdomains","PCASMSetDMSubdomains",osm->dm_subdomains,&osm->dm_subdomains,&flg);CHKERRQ(ierr); 73490d69ab7SBarry Smith ierr = PetscOptionsInt("-pc_asm_blocks","Number of subdomains","PCASMSetTotalSubdomains",osm->n,&blocks,&flg);CHKERRQ(ierr); 73565db9045SDmitry Karpeev if (flg) { 736f77a5242SKarl Rupp ierr = PCASMSetTotalSubdomains(pc,blocks,NULL,NULL);CHKERRQ(ierr); 737d709ea83SDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 73865db9045SDmitry Karpeev } 739342c94f9SMatthew G. Knepley ierr = PetscOptionsInt("-pc_asm_local_blocks","Number of local subdomains","PCASMSetLocalSubdomains",osm->n_local_true,&blocks,&flg);CHKERRQ(ierr); 740342c94f9SMatthew G. Knepley if (flg) { 741342c94f9SMatthew G. Knepley ierr = PCASMSetLocalSubdomains(pc,blocks,NULL,NULL);CHKERRQ(ierr); 742342c94f9SMatthew G. Knepley osm->dm_subdomains = PETSC_FALSE; 743342c94f9SMatthew G. Knepley } 74490d69ab7SBarry Smith ierr = PetscOptionsInt("-pc_asm_overlap","Number of grid points overlap","PCASMSetOverlap",osm->overlap,&ovl,&flg);CHKERRQ(ierr); 74565db9045SDmitry Karpeev if (flg) { 74665db9045SDmitry Karpeev ierr = PCASMSetOverlap(pc,ovl);CHKERRQ(ierr); 747d709ea83SDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 74865db9045SDmitry Karpeev } 74990d69ab7SBarry Smith flg = PETSC_FALSE; 75090d69ab7SBarry Smith ierr = PetscOptionsEnum("-pc_asm_type","Type of restriction/extension","PCASMSetType",PCASMTypes,(PetscEnum)osm->type,(PetscEnum*)&asmtype,&flg);CHKERRQ(ierr); 75192bb6962SLisandro Dalcin if (flg) {ierr = PCASMSetType(pc,asmtype);CHKERRQ(ierr); } 75212cd4985SMatthew G. Knepley flg = PETSC_FALSE; 75312cd4985SMatthew G. Knepley ierr = PetscOptionsEnum("-pc_asm_local_type","Type of local solver composition","PCASMSetLocalType",PCCompositeTypes,(PetscEnum)osm->loctype,(PetscEnum*)&loctype,&flg);CHKERRQ(ierr); 75412cd4985SMatthew G. Knepley if (flg) {ierr = PCASMSetLocalType(pc,loctype);CHKERRQ(ierr); } 75580ec0b7dSPatrick Sanan ierr = PetscOptionsFList("-pc_asm_sub_mat_type","Subsolve Matrix Type","PCASMSetSubMatType",MatList,NULL,sub_mat_type,256,&flg);CHKERRQ(ierr); 75680ec0b7dSPatrick Sanan if (flg) { 757459726d8SSatish Balay ierr = PCASMSetSubMatType(pc,sub_mat_type);CHKERRQ(ierr); 75880ec0b7dSPatrick Sanan } 7594b9ad928SBarry Smith ierr = PetscOptionsTail();CHKERRQ(ierr); 7604b9ad928SBarry Smith PetscFunctionReturn(0); 7614b9ad928SBarry Smith } 7624b9ad928SBarry Smith 7634b9ad928SBarry Smith /*------------------------------------------------------------------------------------*/ 7644b9ad928SBarry Smith 7651e6b0712SBarry Smith static PetscErrorCode PCASMSetLocalSubdomains_ASM(PC pc,PetscInt n,IS is[],IS is_local[]) 7664b9ad928SBarry Smith { 7674b9ad928SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 76892bb6962SLisandro Dalcin PetscErrorCode ierr; 76992bb6962SLisandro Dalcin PetscInt i; 7704b9ad928SBarry Smith 7714b9ad928SBarry Smith PetscFunctionBegin; 772e32f2f54SBarry Smith if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Each process must have 1 or more blocks, n = %D",n); 773ce94432eSBarry 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()."); 774e7e72b3dSBarry Smith 7754b9ad928SBarry Smith if (!pc->setupcalled) { 77692bb6962SLisandro Dalcin if (is) { 77792bb6962SLisandro Dalcin for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)is[i]);CHKERRQ(ierr);} 77892bb6962SLisandro Dalcin } 779832fc9a5SMatthew Knepley if (is_local) { 780832fc9a5SMatthew Knepley for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)is_local[i]);CHKERRQ(ierr);} 781832fc9a5SMatthew Knepley } 7822b691e39SMatthew Knepley ierr = PCASMDestroySubdomains(osm->n_local_true,osm->is,osm->is_local);CHKERRQ(ierr); 7832fa5cd67SKarl Rupp 7844b9ad928SBarry Smith osm->n_local_true = n; 7850a545947SLisandro Dalcin osm->is = NULL; 7860a545947SLisandro Dalcin osm->is_local = NULL; 78792bb6962SLisandro Dalcin if (is) { 788785e854fSJed Brown ierr = PetscMalloc1(n,&osm->is);CHKERRQ(ierr); 7892fa5cd67SKarl Rupp for (i=0; i<n; i++) osm->is[i] = is[i]; 7903d03bb48SJed Brown /* Flag indicating that the user has set overlapping subdomains so PCASM should not increase their size. */ 7913d03bb48SJed Brown osm->overlap = -1; 79292bb6962SLisandro Dalcin } 7932b691e39SMatthew Knepley if (is_local) { 794785e854fSJed Brown ierr = PetscMalloc1(n,&osm->is_local);CHKERRQ(ierr); 7952fa5cd67SKarl Rupp for (i=0; i<n; i++) osm->is_local[i] = is_local[i]; 796a35b7b57SDmitry Karpeev if (!is) { 797785e854fSJed Brown ierr = PetscMalloc1(osm->n_local_true,&osm->is);CHKERRQ(ierr); 798a35b7b57SDmitry Karpeev for (i=0; i<osm->n_local_true; i++) { 799a35b7b57SDmitry Karpeev if (osm->overlap > 0) { /* With positive overlap, osm->is[i] will be modified */ 800a35b7b57SDmitry Karpeev ierr = ISDuplicate(osm->is_local[i],&osm->is[i]);CHKERRQ(ierr); 801a35b7b57SDmitry Karpeev ierr = ISCopy(osm->is_local[i],osm->is[i]);CHKERRQ(ierr); 802a35b7b57SDmitry Karpeev } else { 803a35b7b57SDmitry Karpeev ierr = PetscObjectReference((PetscObject)osm->is_local[i]);CHKERRQ(ierr); 804a35b7b57SDmitry Karpeev osm->is[i] = osm->is_local[i]; 805a35b7b57SDmitry Karpeev } 806a35b7b57SDmitry Karpeev } 807a35b7b57SDmitry Karpeev } 8082b691e39SMatthew Knepley } 8094b9ad928SBarry Smith } 8104b9ad928SBarry Smith PetscFunctionReturn(0); 8114b9ad928SBarry Smith } 8124b9ad928SBarry Smith 8131e6b0712SBarry Smith static PetscErrorCode PCASMSetTotalSubdomains_ASM(PC pc,PetscInt N,IS *is,IS *is_local) 8144b9ad928SBarry Smith { 8154b9ad928SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 8166849ba73SBarry Smith PetscErrorCode ierr; 81713f74950SBarry Smith PetscMPIInt rank,size; 81878904715SLisandro Dalcin PetscInt n; 8194b9ad928SBarry Smith 8204b9ad928SBarry Smith PetscFunctionBegin; 821ce94432eSBarry Smith if (N < 1) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Number of total blocks must be > 0, N = %D",N); 822ce94432eSBarry 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."); 8234b9ad928SBarry Smith 8244b9ad928SBarry Smith /* 825880770e9SJed Brown Split the subdomains equally among all processors 8264b9ad928SBarry Smith */ 827ce94432eSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);CHKERRQ(ierr); 828ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr); 8294b9ad928SBarry Smith n = N/size + ((N % size) > rank); 830e32f2f54SBarry 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); 831e32f2f54SBarry Smith if (pc->setupcalled && n != osm->n_local_true) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"PCASMSetTotalSubdomains() should be called before PCSetUp()."); 8324b9ad928SBarry Smith if (!pc->setupcalled) { 8332b691e39SMatthew Knepley ierr = PCASMDestroySubdomains(osm->n_local_true,osm->is,osm->is_local);CHKERRQ(ierr); 8342fa5cd67SKarl Rupp 8354b9ad928SBarry Smith osm->n_local_true = n; 8360a545947SLisandro Dalcin osm->is = NULL; 8370a545947SLisandro Dalcin osm->is_local = NULL; 8384b9ad928SBarry Smith } 8394b9ad928SBarry Smith PetscFunctionReturn(0); 8404b9ad928SBarry Smith } 8414b9ad928SBarry Smith 8421e6b0712SBarry Smith static PetscErrorCode PCASMSetOverlap_ASM(PC pc,PetscInt ovl) 8434b9ad928SBarry Smith { 84492bb6962SLisandro Dalcin PC_ASM *osm = (PC_ASM*)pc->data; 8454b9ad928SBarry Smith 8464b9ad928SBarry Smith PetscFunctionBegin; 847ce94432eSBarry Smith if (ovl < 0) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap value requested"); 848ce94432eSBarry Smith if (pc->setupcalled && ovl != osm->overlap) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCASMSetOverlap() should be called before PCSetUp()."); 8492fa5cd67SKarl Rupp if (!pc->setupcalled) osm->overlap = ovl; 8504b9ad928SBarry Smith PetscFunctionReturn(0); 8514b9ad928SBarry Smith } 8524b9ad928SBarry Smith 8531e6b0712SBarry Smith static PetscErrorCode PCASMSetType_ASM(PC pc,PCASMType type) 8544b9ad928SBarry Smith { 85592bb6962SLisandro Dalcin PC_ASM *osm = (PC_ASM*)pc->data; 8564b9ad928SBarry Smith 8574b9ad928SBarry Smith PetscFunctionBegin; 8584b9ad928SBarry Smith osm->type = type; 859bf108f30SBarry Smith osm->type_set = PETSC_TRUE; 8604b9ad928SBarry Smith PetscFunctionReturn(0); 8614b9ad928SBarry Smith } 8624b9ad928SBarry Smith 863c60c7ad4SBarry Smith static PetscErrorCode PCASMGetType_ASM(PC pc,PCASMType *type) 864c60c7ad4SBarry Smith { 865c60c7ad4SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 866c60c7ad4SBarry Smith 867c60c7ad4SBarry Smith PetscFunctionBegin; 868c60c7ad4SBarry Smith *type = osm->type; 869c60c7ad4SBarry Smith PetscFunctionReturn(0); 870c60c7ad4SBarry Smith } 871c60c7ad4SBarry Smith 87212cd4985SMatthew G. Knepley static PetscErrorCode PCASMSetLocalType_ASM(PC pc, PCCompositeType type) 87312cd4985SMatthew G. Knepley { 87412cd4985SMatthew G. Knepley PC_ASM *osm = (PC_ASM *) pc->data; 87512cd4985SMatthew G. Knepley 87612cd4985SMatthew G. Knepley PetscFunctionBegin; 877d0ecd4dfSBarry 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"); 87812cd4985SMatthew G. Knepley osm->loctype = type; 87912cd4985SMatthew G. Knepley PetscFunctionReturn(0); 88012cd4985SMatthew G. Knepley } 88112cd4985SMatthew G. Knepley 88212cd4985SMatthew G. Knepley static PetscErrorCode PCASMGetLocalType_ASM(PC pc, PCCompositeType *type) 88312cd4985SMatthew G. Knepley { 88412cd4985SMatthew G. Knepley PC_ASM *osm = (PC_ASM *) pc->data; 88512cd4985SMatthew G. Knepley 88612cd4985SMatthew G. Knepley PetscFunctionBegin; 88712cd4985SMatthew G. Knepley *type = osm->loctype; 88812cd4985SMatthew G. Knepley PetscFunctionReturn(0); 88912cd4985SMatthew G. Knepley } 89012cd4985SMatthew G. Knepley 8911e6b0712SBarry Smith static PetscErrorCode PCASMSetSortIndices_ASM(PC pc,PetscBool doSort) 8926ed231c7SMatthew Knepley { 8936ed231c7SMatthew Knepley PC_ASM *osm = (PC_ASM*)pc->data; 8946ed231c7SMatthew Knepley 8956ed231c7SMatthew Knepley PetscFunctionBegin; 8966ed231c7SMatthew Knepley osm->sort_indices = doSort; 8976ed231c7SMatthew Knepley PetscFunctionReturn(0); 8986ed231c7SMatthew Knepley } 8996ed231c7SMatthew Knepley 9001e6b0712SBarry Smith static PetscErrorCode PCASMGetSubKSP_ASM(PC pc,PetscInt *n_local,PetscInt *first_local,KSP **ksp) 9014b9ad928SBarry Smith { 90292bb6962SLisandro Dalcin PC_ASM *osm = (PC_ASM*)pc->data; 903dfbe8321SBarry Smith PetscErrorCode ierr; 9044b9ad928SBarry Smith 9054b9ad928SBarry Smith PetscFunctionBegin; 90634a84908Sprj- 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"); 9074b9ad928SBarry Smith 9082fa5cd67SKarl Rupp if (n_local) *n_local = osm->n_local_true; 90992bb6962SLisandro Dalcin if (first_local) { 910ce94432eSBarry Smith ierr = MPI_Scan(&osm->n_local_true,first_local,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr); 91192bb6962SLisandro Dalcin *first_local -= osm->n_local_true; 91292bb6962SLisandro Dalcin } 91392bb6962SLisandro Dalcin if (ksp) { 91492bb6962SLisandro Dalcin /* Assume that local solves are now different; not necessarily 91592bb6962SLisandro Dalcin true though! This flag is used only for PCView_ASM() */ 91692bb6962SLisandro Dalcin *ksp = osm->ksp; 91792bb6962SLisandro Dalcin osm->same_local_solves = PETSC_FALSE; 91892bb6962SLisandro Dalcin } 9194b9ad928SBarry Smith PetscFunctionReturn(0); 9204b9ad928SBarry Smith } 9214b9ad928SBarry Smith 92280ec0b7dSPatrick Sanan static PetscErrorCode PCASMGetSubMatType_ASM(PC pc,MatType *sub_mat_type) 92380ec0b7dSPatrick Sanan { 92480ec0b7dSPatrick Sanan PC_ASM *osm = (PC_ASM*)pc->data; 92580ec0b7dSPatrick Sanan 92680ec0b7dSPatrick Sanan PetscFunctionBegin; 92780ec0b7dSPatrick Sanan PetscValidHeaderSpecific(pc,PC_CLASSID,1); 92880ec0b7dSPatrick Sanan PetscValidPointer(sub_mat_type,2); 92980ec0b7dSPatrick Sanan *sub_mat_type = osm->sub_mat_type; 93080ec0b7dSPatrick Sanan PetscFunctionReturn(0); 93180ec0b7dSPatrick Sanan } 93280ec0b7dSPatrick Sanan 93380ec0b7dSPatrick Sanan static PetscErrorCode PCASMSetSubMatType_ASM(PC pc,MatType sub_mat_type) 93480ec0b7dSPatrick Sanan { 93580ec0b7dSPatrick Sanan PetscErrorCode ierr; 93680ec0b7dSPatrick Sanan PC_ASM *osm = (PC_ASM*)pc->data; 93780ec0b7dSPatrick Sanan 93880ec0b7dSPatrick Sanan PetscFunctionBegin; 93980ec0b7dSPatrick Sanan PetscValidHeaderSpecific(pc,PC_CLASSID,1); 94080ec0b7dSPatrick Sanan ierr = PetscFree(osm->sub_mat_type);CHKERRQ(ierr); 94180ec0b7dSPatrick Sanan ierr = PetscStrallocpy(sub_mat_type,(char**)&osm->sub_mat_type);CHKERRQ(ierr); 94280ec0b7dSPatrick Sanan PetscFunctionReturn(0); 94380ec0b7dSPatrick Sanan } 94480ec0b7dSPatrick Sanan 9454b9ad928SBarry Smith /*@C 9461093a601SBarry Smith PCASMSetLocalSubdomains - Sets the local subdomains (for this processor only) for the additive Schwarz preconditioner. 9474b9ad928SBarry Smith 948d083f849SBarry Smith Collective on pc 9494b9ad928SBarry Smith 9504b9ad928SBarry Smith Input Parameters: 9514b9ad928SBarry Smith + pc - the preconditioner context 9524b9ad928SBarry Smith . n - the number of subdomains for this processor (default value = 1) 9538c03b21aSDmitry Karpeev . is - the index set that defines the subdomains for this processor 9540298fd71SBarry Smith (or NULL for PETSc to determine subdomains) 955f1ee410cSBarry 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 956f1ee410cSBarry Smith (or NULL to not provide these) 9574b9ad928SBarry Smith 958342c94f9SMatthew G. Knepley Options Database Key: 959342c94f9SMatthew G. Knepley To set the total number of subdomain blocks rather than specify the 960342c94f9SMatthew G. Knepley index sets, use the option 961342c94f9SMatthew G. Knepley . -pc_asm_local_blocks <blks> - Sets local blocks 962342c94f9SMatthew G. Knepley 9634b9ad928SBarry Smith Notes: 9641093a601SBarry Smith The IS numbering is in the parallel, global numbering of the vector for both is and is_local 9654b9ad928SBarry Smith 9664b9ad928SBarry Smith By default the ASM preconditioner uses 1 block per processor. 9674b9ad928SBarry Smith 9684b9ad928SBarry Smith Use PCASMSetTotalSubdomains() to set the subdomains for all processors. 9694b9ad928SBarry Smith 970f1ee410cSBarry Smith If is_local is provided and PCASMType is PC_ASM_RESTRICT then the solution only over the is_local region is interpolated 971f1ee410cSBarry Smith back to form the global solution (this is the standard restricted additive Schwarz method) 972f1ee410cSBarry Smith 973f1ee410cSBarry 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 974f1ee410cSBarry Smith no code to handle that case. 975f1ee410cSBarry Smith 9764b9ad928SBarry Smith Level: advanced 9774b9ad928SBarry Smith 9784b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(), 979f1ee410cSBarry Smith PCASMCreateSubdomains2D(), PCASMGetLocalSubdomains(), PCASMType, PCASMSetType() 9804b9ad928SBarry Smith @*/ 9817087cfbeSBarry Smith PetscErrorCode PCASMSetLocalSubdomains(PC pc,PetscInt n,IS is[],IS is_local[]) 9824b9ad928SBarry Smith { 9837bb14e67SBarry Smith PetscErrorCode ierr; 9844b9ad928SBarry Smith 9854b9ad928SBarry Smith PetscFunctionBegin; 9860700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 9877bb14e67SBarry Smith ierr = PetscTryMethod(pc,"PCASMSetLocalSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,n,is,is_local));CHKERRQ(ierr); 9884b9ad928SBarry Smith PetscFunctionReturn(0); 9894b9ad928SBarry Smith } 9904b9ad928SBarry Smith 9914b9ad928SBarry Smith /*@C 992feb221f8SDmitry Karpeev PCASMSetTotalSubdomains - Sets the subdomains for all processors for the 9934b9ad928SBarry Smith additive Schwarz preconditioner. Either all or no processors in the 9944b9ad928SBarry Smith PC communicator must call this routine, with the same index sets. 9954b9ad928SBarry Smith 996d083f849SBarry Smith Collective on pc 9974b9ad928SBarry Smith 9984b9ad928SBarry Smith Input Parameters: 9994b9ad928SBarry Smith + pc - the preconditioner context 1000feb221f8SDmitry Karpeev . N - the number of subdomains for all processors 1001feb221f8SDmitry Karpeev . is - the index sets that define the subdomains for all processors 1002dfaa0c5fSBarry Smith (or NULL to ask PETSc to determine the subdomains) 10032b691e39SMatthew Knepley - is_local - the index sets that define the local part of the subdomains for this processor 1004f1ee410cSBarry Smith (or NULL to not provide this information) 10054b9ad928SBarry Smith 10064b9ad928SBarry Smith Options Database Key: 10074b9ad928SBarry Smith To set the total number of subdomain blocks rather than specify the 10084b9ad928SBarry Smith index sets, use the option 10094b9ad928SBarry Smith . -pc_asm_blocks <blks> - Sets total blocks 10104b9ad928SBarry Smith 10114b9ad928SBarry Smith Notes: 1012f1ee410cSBarry Smith Currently you cannot use this to set the actual subdomains with the argument is or is_local. 10134b9ad928SBarry Smith 10144b9ad928SBarry Smith By default the ASM preconditioner uses 1 block per processor. 10154b9ad928SBarry Smith 10164b9ad928SBarry Smith These index sets cannot be destroyed until after completion of the 10174b9ad928SBarry Smith linear solves for which the ASM preconditioner is being used. 10184b9ad928SBarry Smith 10194b9ad928SBarry Smith Use PCASMSetLocalSubdomains() to set local subdomains. 10204b9ad928SBarry Smith 10211093a601SBarry Smith The IS numbering is in the parallel, global numbering of the vector for both is and is_local 10221093a601SBarry Smith 10234b9ad928SBarry Smith Level: advanced 10244b9ad928SBarry Smith 10254b9ad928SBarry Smith .seealso: PCASMSetLocalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(), 10264b9ad928SBarry Smith PCASMCreateSubdomains2D() 10274b9ad928SBarry Smith @*/ 10287087cfbeSBarry Smith PetscErrorCode PCASMSetTotalSubdomains(PC pc,PetscInt N,IS is[],IS is_local[]) 10294b9ad928SBarry Smith { 10307bb14e67SBarry Smith PetscErrorCode ierr; 10314b9ad928SBarry Smith 10324b9ad928SBarry Smith PetscFunctionBegin; 10330700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 10347bb14e67SBarry Smith ierr = PetscTryMethod(pc,"PCASMSetTotalSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,N,is,is_local));CHKERRQ(ierr); 10354b9ad928SBarry Smith PetscFunctionReturn(0); 10364b9ad928SBarry Smith } 10374b9ad928SBarry Smith 10384b9ad928SBarry Smith /*@ 10394b9ad928SBarry Smith PCASMSetOverlap - Sets the overlap between a pair of subdomains for the 10404b9ad928SBarry Smith additive Schwarz preconditioner. Either all or no processors in the 1041f1ee410cSBarry Smith PC communicator must call this routine. 10424b9ad928SBarry Smith 1043d083f849SBarry Smith Logically Collective on pc 10444b9ad928SBarry Smith 10454b9ad928SBarry Smith Input Parameters: 10464b9ad928SBarry Smith + pc - the preconditioner context 10474b9ad928SBarry Smith - ovl - the amount of overlap between subdomains (ovl >= 0, default value = 1) 10484b9ad928SBarry Smith 10494b9ad928SBarry Smith Options Database Key: 10504b9ad928SBarry Smith . -pc_asm_overlap <ovl> - Sets overlap 10514b9ad928SBarry Smith 10524b9ad928SBarry Smith Notes: 10534b9ad928SBarry Smith By default the ASM preconditioner uses 1 block per processor. To use 10544b9ad928SBarry Smith multiple blocks per perocessor, see PCASMSetTotalSubdomains() and 10554b9ad928SBarry Smith PCASMSetLocalSubdomains() (and the option -pc_asm_blocks <blks>). 10564b9ad928SBarry Smith 10574b9ad928SBarry Smith The overlap defaults to 1, so if one desires that no additional 10584b9ad928SBarry Smith overlap be computed beyond what may have been set with a call to 10594b9ad928SBarry Smith PCASMSetTotalSubdomains() or PCASMSetLocalSubdomains(), then ovl 10604b9ad928SBarry Smith must be set to be 0. In particular, if one does not explicitly set 10614b9ad928SBarry Smith the subdomains an application code, then all overlap would be computed 10624b9ad928SBarry Smith internally by PETSc, and using an overlap of 0 would result in an ASM 10634b9ad928SBarry Smith variant that is equivalent to the block Jacobi preconditioner. 10644b9ad928SBarry Smith 1065f1ee410cSBarry Smith The default algorithm used by PETSc to increase overlap is fast, but not scalable, 1066f1ee410cSBarry Smith use the option -mat_increase_overlap_scalable when the problem and number of processes is large. 1067f1ee410cSBarry Smith 10684b9ad928SBarry Smith Note that one can define initial index sets with any overlap via 1069f1ee410cSBarry Smith PCASMSetLocalSubdomains(); the routine 10704b9ad928SBarry Smith PCASMSetOverlap() merely allows PETSc to extend that overlap further 10714b9ad928SBarry Smith if desired. 10724b9ad928SBarry Smith 10734b9ad928SBarry Smith Level: intermediate 10744b9ad928SBarry Smith 10754b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetLocalSubdomains(), PCASMGetSubKSP(), 1076f1ee410cSBarry Smith PCASMCreateSubdomains2D(), PCASMGetLocalSubdomains(), MatIncreaseOverlap() 10774b9ad928SBarry Smith @*/ 10787087cfbeSBarry Smith PetscErrorCode PCASMSetOverlap(PC pc,PetscInt ovl) 10794b9ad928SBarry Smith { 10807bb14e67SBarry Smith PetscErrorCode ierr; 10814b9ad928SBarry Smith 10824b9ad928SBarry Smith PetscFunctionBegin; 10830700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1084c5eb9154SBarry Smith PetscValidLogicalCollectiveInt(pc,ovl,2); 10857bb14e67SBarry Smith ierr = PetscTryMethod(pc,"PCASMSetOverlap_C",(PC,PetscInt),(pc,ovl));CHKERRQ(ierr); 10864b9ad928SBarry Smith PetscFunctionReturn(0); 10874b9ad928SBarry Smith } 10884b9ad928SBarry Smith 10894b9ad928SBarry Smith /*@ 10904b9ad928SBarry Smith PCASMSetType - Sets the type of restriction and interpolation used 10914b9ad928SBarry Smith for local problems in the additive Schwarz method. 10924b9ad928SBarry Smith 1093d083f849SBarry Smith Logically Collective on pc 10944b9ad928SBarry Smith 10954b9ad928SBarry Smith Input Parameters: 10964b9ad928SBarry Smith + pc - the preconditioner context 10974b9ad928SBarry Smith - type - variant of ASM, one of 10984b9ad928SBarry Smith .vb 10994b9ad928SBarry Smith PC_ASM_BASIC - full interpolation and restriction 11004b9ad928SBarry Smith PC_ASM_RESTRICT - full restriction, local processor interpolation 11014b9ad928SBarry Smith PC_ASM_INTERPOLATE - full interpolation, local processor restriction 11024b9ad928SBarry Smith PC_ASM_NONE - local processor restriction and interpolation 11034b9ad928SBarry Smith .ve 11044b9ad928SBarry Smith 11054b9ad928SBarry Smith Options Database Key: 11064b9ad928SBarry Smith . -pc_asm_type [basic,restrict,interpolate,none] - Sets ASM type 11074b9ad928SBarry Smith 110895452b02SPatrick Sanan Notes: 110995452b02SPatrick Sanan if the is_local arguments are passed to PCASMSetLocalSubdomains() then they are used when PC_ASM_RESTRICT has been selected 1110f1ee410cSBarry Smith to limit the local processor interpolation 1111f1ee410cSBarry Smith 11124b9ad928SBarry Smith Level: intermediate 11134b9ad928SBarry Smith 11144b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetTotalSubdomains(), PCASMGetSubKSP(), 1115f1ee410cSBarry Smith PCASMCreateSubdomains2D(), PCASMType, PCASMSetLocalType(), PCASMGetLocalType() 11164b9ad928SBarry Smith @*/ 11177087cfbeSBarry Smith PetscErrorCode PCASMSetType(PC pc,PCASMType type) 11184b9ad928SBarry Smith { 11197bb14e67SBarry Smith PetscErrorCode ierr; 11204b9ad928SBarry Smith 11214b9ad928SBarry Smith PetscFunctionBegin; 11220700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1123c5eb9154SBarry Smith PetscValidLogicalCollectiveEnum(pc,type,2); 11247bb14e67SBarry Smith ierr = PetscTryMethod(pc,"PCASMSetType_C",(PC,PCASMType),(pc,type));CHKERRQ(ierr); 11254b9ad928SBarry Smith PetscFunctionReturn(0); 11264b9ad928SBarry Smith } 11274b9ad928SBarry Smith 1128c60c7ad4SBarry Smith /*@ 1129c60c7ad4SBarry Smith PCASMGetType - Gets the type of restriction and interpolation used 1130c60c7ad4SBarry Smith for local problems in the additive Schwarz method. 1131c60c7ad4SBarry Smith 1132d083f849SBarry Smith Logically Collective on pc 1133c60c7ad4SBarry Smith 1134c60c7ad4SBarry Smith Input Parameter: 1135c60c7ad4SBarry Smith . pc - the preconditioner context 1136c60c7ad4SBarry Smith 1137c60c7ad4SBarry Smith Output Parameter: 1138c60c7ad4SBarry Smith . type - variant of ASM, one of 1139c60c7ad4SBarry Smith 1140c60c7ad4SBarry Smith .vb 1141c60c7ad4SBarry Smith PC_ASM_BASIC - full interpolation and restriction 1142c60c7ad4SBarry Smith PC_ASM_RESTRICT - full restriction, local processor interpolation 1143c60c7ad4SBarry Smith PC_ASM_INTERPOLATE - full interpolation, local processor restriction 1144c60c7ad4SBarry Smith PC_ASM_NONE - local processor restriction and interpolation 1145c60c7ad4SBarry Smith .ve 1146c60c7ad4SBarry Smith 1147c60c7ad4SBarry Smith Options Database Key: 1148c60c7ad4SBarry Smith . -pc_asm_type [basic,restrict,interpolate,none] - Sets ASM type 1149c60c7ad4SBarry Smith 1150c60c7ad4SBarry Smith Level: intermediate 1151c60c7ad4SBarry Smith 1152c60c7ad4SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetTotalSubdomains(), PCASMGetSubKSP(), 1153f1ee410cSBarry Smith PCASMCreateSubdomains2D(), PCASMType, PCASMSetType(), PCASMSetLocalType(), PCASMGetLocalType() 1154c60c7ad4SBarry Smith @*/ 1155c60c7ad4SBarry Smith PetscErrorCode PCASMGetType(PC pc,PCASMType *type) 1156c60c7ad4SBarry Smith { 1157c60c7ad4SBarry Smith PetscErrorCode ierr; 1158c60c7ad4SBarry Smith 1159c60c7ad4SBarry Smith PetscFunctionBegin; 1160c60c7ad4SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1161c60c7ad4SBarry Smith ierr = PetscUseMethod(pc,"PCASMGetType_C",(PC,PCASMType*),(pc,type));CHKERRQ(ierr); 1162c60c7ad4SBarry Smith PetscFunctionReturn(0); 1163c60c7ad4SBarry Smith } 1164c60c7ad4SBarry Smith 116512cd4985SMatthew G. Knepley /*@ 116612cd4985SMatthew G. Knepley PCASMSetLocalType - Sets the type of composition used for local problems in the additive Schwarz method. 116712cd4985SMatthew G. Knepley 1168d083f849SBarry Smith Logically Collective on pc 116912cd4985SMatthew G. Knepley 117012cd4985SMatthew G. Knepley Input Parameters: 117112cd4985SMatthew G. Knepley + pc - the preconditioner context 117212cd4985SMatthew G. Knepley - type - type of composition, one of 117312cd4985SMatthew G. Knepley .vb 117412cd4985SMatthew G. Knepley PC_COMPOSITE_ADDITIVE - local additive combination 117512cd4985SMatthew G. Knepley PC_COMPOSITE_MULTIPLICATIVE - local multiplicative combination 117612cd4985SMatthew G. Knepley .ve 117712cd4985SMatthew G. Knepley 117812cd4985SMatthew G. Knepley Options Database Key: 117912cd4985SMatthew G. Knepley . -pc_asm_local_type [additive,multiplicative] - Sets local solver composition type 118012cd4985SMatthew G. Knepley 118112cd4985SMatthew G. Knepley Level: intermediate 118212cd4985SMatthew G. Knepley 1183f1ee410cSBarry Smith .seealso: PCASMSetType(), PCASMGetType(), PCASMGetLocalType(), PCASM, PCASMType, PCASMSetType(), PCASMGetType(), PCCompositeType 118412cd4985SMatthew G. Knepley @*/ 118512cd4985SMatthew G. Knepley PetscErrorCode PCASMSetLocalType(PC pc, PCCompositeType type) 118612cd4985SMatthew G. Knepley { 118712cd4985SMatthew G. Knepley PetscErrorCode ierr; 118812cd4985SMatthew G. Knepley 118912cd4985SMatthew G. Knepley PetscFunctionBegin; 119012cd4985SMatthew G. Knepley PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 119112cd4985SMatthew G. Knepley PetscValidLogicalCollectiveEnum(pc, type, 2); 119212cd4985SMatthew G. Knepley ierr = PetscTryMethod(pc, "PCASMSetLocalType_C", (PC, PCCompositeType), (pc, type));CHKERRQ(ierr); 119312cd4985SMatthew G. Knepley PetscFunctionReturn(0); 119412cd4985SMatthew G. Knepley } 119512cd4985SMatthew G. Knepley 119612cd4985SMatthew G. Knepley /*@ 119712cd4985SMatthew G. Knepley PCASMGetLocalType - Gets the type of composition used for local problems in the additive Schwarz method. 119812cd4985SMatthew G. Knepley 1199d083f849SBarry Smith Logically Collective on pc 120012cd4985SMatthew G. Knepley 120112cd4985SMatthew G. Knepley Input Parameter: 120212cd4985SMatthew G. Knepley . pc - the preconditioner context 120312cd4985SMatthew G. Knepley 120412cd4985SMatthew G. Knepley Output Parameter: 120512cd4985SMatthew G. Knepley . type - type of composition, one of 120612cd4985SMatthew G. Knepley .vb 120712cd4985SMatthew G. Knepley PC_COMPOSITE_ADDITIVE - local additive combination 120812cd4985SMatthew G. Knepley PC_COMPOSITE_MULTIPLICATIVE - local multiplicative combination 120912cd4985SMatthew G. Knepley .ve 121012cd4985SMatthew G. Knepley 121112cd4985SMatthew G. Knepley Options Database Key: 121212cd4985SMatthew G. Knepley . -pc_asm_local_type [additive,multiplicative] - Sets local solver composition type 121312cd4985SMatthew G. Knepley 121412cd4985SMatthew G. Knepley Level: intermediate 121512cd4985SMatthew G. Knepley 1216f1ee410cSBarry Smith .seealso: PCASMSetType(), PCASMGetType(), PCASMSetLocalType(), PCASMCreate(), PCASMType, PCASMSetType(), PCASMGetType(), PCCompositeType 121712cd4985SMatthew G. Knepley @*/ 121812cd4985SMatthew G. Knepley PetscErrorCode PCASMGetLocalType(PC pc, PCCompositeType *type) 121912cd4985SMatthew G. Knepley { 122012cd4985SMatthew G. Knepley PetscErrorCode ierr; 122112cd4985SMatthew G. Knepley 122212cd4985SMatthew G. Knepley PetscFunctionBegin; 122312cd4985SMatthew G. Knepley PetscValidHeaderSpecific(pc, PC_CLASSID, 1); 122412cd4985SMatthew G. Knepley PetscValidPointer(type, 2); 122512cd4985SMatthew G. Knepley ierr = PetscUseMethod(pc, "PCASMGetLocalType_C", (PC, PCCompositeType *), (pc, type));CHKERRQ(ierr); 122612cd4985SMatthew G. Knepley PetscFunctionReturn(0); 122712cd4985SMatthew G. Knepley } 122812cd4985SMatthew G. Knepley 12296ed231c7SMatthew Knepley /*@ 12306ed231c7SMatthew Knepley PCASMSetSortIndices - Determines whether subdomain indices are sorted. 12316ed231c7SMatthew Knepley 1232d083f849SBarry Smith Logically Collective on pc 12336ed231c7SMatthew Knepley 12346ed231c7SMatthew Knepley Input Parameters: 12356ed231c7SMatthew Knepley + pc - the preconditioner context 12366ed231c7SMatthew Knepley - doSort - sort the subdomain indices 12376ed231c7SMatthew Knepley 12386ed231c7SMatthew Knepley Level: intermediate 12396ed231c7SMatthew Knepley 12406ed231c7SMatthew Knepley .seealso: PCASMSetLocalSubdomains(), PCASMSetTotalSubdomains(), PCASMGetSubKSP(), 12416ed231c7SMatthew Knepley PCASMCreateSubdomains2D() 12426ed231c7SMatthew Knepley @*/ 12437087cfbeSBarry Smith PetscErrorCode PCASMSetSortIndices(PC pc,PetscBool doSort) 12446ed231c7SMatthew Knepley { 12457bb14e67SBarry Smith PetscErrorCode ierr; 12466ed231c7SMatthew Knepley 12476ed231c7SMatthew Knepley PetscFunctionBegin; 12480700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1249acfcf0e5SJed Brown PetscValidLogicalCollectiveBool(pc,doSort,2); 12507bb14e67SBarry Smith ierr = PetscTryMethod(pc,"PCASMSetSortIndices_C",(PC,PetscBool),(pc,doSort));CHKERRQ(ierr); 12516ed231c7SMatthew Knepley PetscFunctionReturn(0); 12526ed231c7SMatthew Knepley } 12536ed231c7SMatthew Knepley 12544b9ad928SBarry Smith /*@C 12554b9ad928SBarry Smith PCASMGetSubKSP - Gets the local KSP contexts for all blocks on 12564b9ad928SBarry Smith this processor. 12574b9ad928SBarry Smith 1258d083f849SBarry Smith Collective on pc iff first_local is requested 12594b9ad928SBarry Smith 12604b9ad928SBarry Smith Input Parameter: 12614b9ad928SBarry Smith . pc - the preconditioner context 12624b9ad928SBarry Smith 12634b9ad928SBarry Smith Output Parameters: 12640298fd71SBarry Smith + n_local - the number of blocks on this processor or NULL 12650298fd71SBarry Smith . first_local - the global number of the first block on this processor or NULL, 12660298fd71SBarry Smith all processors must request or all must pass NULL 12674b9ad928SBarry Smith - ksp - the array of KSP contexts 12684b9ad928SBarry Smith 12694b9ad928SBarry Smith Note: 1270d29017ddSJed Brown After PCASMGetSubKSP() the array of KSPes is not to be freed. 12714b9ad928SBarry Smith 12724b9ad928SBarry Smith You must call KSPSetUp() before calling PCASMGetSubKSP(). 12734b9ad928SBarry Smith 1274d29017ddSJed Brown Fortran note: 12752bf68e3eSBarry 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. 1276d29017ddSJed Brown 12774b9ad928SBarry Smith Level: advanced 12784b9ad928SBarry Smith 12794b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetTotalSubdomains(), PCASMSetOverlap(), 12804b9ad928SBarry Smith PCASMCreateSubdomains2D(), 12814b9ad928SBarry Smith @*/ 12827087cfbeSBarry Smith PetscErrorCode PCASMGetSubKSP(PC pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp[]) 12834b9ad928SBarry Smith { 12847bb14e67SBarry Smith PetscErrorCode ierr; 12854b9ad928SBarry Smith 12864b9ad928SBarry Smith PetscFunctionBegin; 12870700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 12887bb14e67SBarry Smith ierr = PetscUseMethod(pc,"PCASMGetSubKSP_C",(PC,PetscInt*,PetscInt*,KSP **),(pc,n_local,first_local,ksp));CHKERRQ(ierr); 12894b9ad928SBarry Smith PetscFunctionReturn(0); 12904b9ad928SBarry Smith } 12914b9ad928SBarry Smith 12924b9ad928SBarry Smith /* -------------------------------------------------------------------------------------*/ 12934b9ad928SBarry Smith /*MC 12943b09bd56SBarry Smith PCASM - Use the (restricted) additive Schwarz method, each block is (approximately) solved with 12954b9ad928SBarry Smith its own KSP object. 12964b9ad928SBarry Smith 12974b9ad928SBarry Smith Options Database Keys: 129849517cdeSBarry Smith + -pc_asm_blocks <blks> - Sets total blocks 12994b9ad928SBarry Smith . -pc_asm_overlap <ovl> - Sets overlap 1300f1ee410cSBarry Smith . -pc_asm_type [basic,restrict,interpolate,none] - Sets ASM type, default is restrict 1301f1ee410cSBarry Smith - -pc_asm_local_type [additive, multiplicative] - Sets ASM type, default is additive 13024b9ad928SBarry Smith 13033b09bd56SBarry Smith IMPORTANT: If you run with, for example, 3 blocks on 1 processor or 3 blocks on 3 processors you 13043b09bd56SBarry Smith will get a different convergence rate due to the default option of -pc_asm_type restrict. Use 13053b09bd56SBarry Smith -pc_asm_type basic to use the standard ASM. 13063b09bd56SBarry Smith 130795452b02SPatrick Sanan Notes: 130895452b02SPatrick Sanan Each processor can have one or more blocks, but a block cannot be shared by more 1309f1ee410cSBarry Smith than one processor. Use PCGASM for subdomains shared by multiple processes. Defaults to one block per processor. 13104b9ad928SBarry Smith 13113b09bd56SBarry Smith To set options on the solvers for each block append -sub_ to all the KSP, and PC 1312d7ee0231SBarry Smith options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly 13134b9ad928SBarry Smith 1314a8c7a070SBarry Smith To set the options on the solvers separate for each block call PCASMGetSubKSP() 13154b9ad928SBarry Smith and set the options directly on the resulting KSP object (you can access its PC 13164b9ad928SBarry Smith with KSPGetPC()) 13174b9ad928SBarry Smith 13184b9ad928SBarry Smith Level: beginner 13194b9ad928SBarry Smith 1320c582cd25SBarry Smith References: 132196a0c994SBarry Smith + 1. - M Dryja, OB Widlund, An additive variant of the Schwarz alternating method for the case of many subregions 132296a0c994SBarry Smith Courant Institute, New York University Technical report 13236d33885cSprj- - 2. - Barry Smith, Petter Bjorstad, and William Gropp, Domain Decompositions: Parallel Multilevel Methods for Elliptic Partial Differential Equations, 132496a0c994SBarry Smith Cambridge University Press. 1325c582cd25SBarry Smith 13264b9ad928SBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 1327f1ee410cSBarry Smith PCBJACOBI, PCASMGetSubKSP(), PCASMSetLocalSubdomains(), PCASMType, PCASMGetType(), PCASMSetLocalType(), PCASMGetLocalType() 132834a84908Sprj- PCASMSetTotalSubdomains(), PCSetModifySubMatrices(), PCASMSetOverlap(), PCASMSetType(), PCCompositeType 1329e09e08ccSBarry Smith 13304b9ad928SBarry Smith M*/ 13314b9ad928SBarry Smith 13328cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_ASM(PC pc) 13334b9ad928SBarry Smith { 1334dfbe8321SBarry Smith PetscErrorCode ierr; 13354b9ad928SBarry Smith PC_ASM *osm; 13364b9ad928SBarry Smith 13374b9ad928SBarry Smith PetscFunctionBegin; 1338b00a9115SJed Brown ierr = PetscNewLog(pc,&osm);CHKERRQ(ierr); 13392fa5cd67SKarl Rupp 13404b9ad928SBarry Smith osm->n = PETSC_DECIDE; 13414b9ad928SBarry Smith osm->n_local = 0; 13422b837212SDmitry Karpeev osm->n_local_true = PETSC_DECIDE; 13434b9ad928SBarry Smith osm->overlap = 1; 13440a545947SLisandro Dalcin osm->ksp = NULL; 13450a545947SLisandro Dalcin osm->restriction = NULL; 13460a545947SLisandro Dalcin osm->lprolongation = NULL; 13470a545947SLisandro Dalcin osm->lrestriction = NULL; 13480a545947SLisandro Dalcin osm->x = NULL; 13490a545947SLisandro Dalcin osm->y = NULL; 13500a545947SLisandro Dalcin osm->is = NULL; 13510a545947SLisandro Dalcin osm->is_local = NULL; 13520a545947SLisandro Dalcin osm->mat = NULL; 13530a545947SLisandro Dalcin osm->pmat = NULL; 13544b9ad928SBarry Smith osm->type = PC_ASM_RESTRICT; 135512cd4985SMatthew G. Knepley osm->loctype = PC_COMPOSITE_ADDITIVE; 13564b9ad928SBarry Smith osm->same_local_solves = PETSC_TRUE; 13576ed231c7SMatthew Knepley osm->sort_indices = PETSC_TRUE; 1358d709ea83SDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 135980ec0b7dSPatrick Sanan osm->sub_mat_type = NULL; 13604b9ad928SBarry Smith 136192bb6962SLisandro Dalcin pc->data = (void*)osm; 13624b9ad928SBarry Smith pc->ops->apply = PCApply_ASM; 136348e38a8aSPierre Jolivet pc->ops->matapply = PCMatApply_ASM; 13644b9ad928SBarry Smith pc->ops->applytranspose = PCApplyTranspose_ASM; 13654b9ad928SBarry Smith pc->ops->setup = PCSetUp_ASM; 1366e91c6855SBarry Smith pc->ops->reset = PCReset_ASM; 13674b9ad928SBarry Smith pc->ops->destroy = PCDestroy_ASM; 13684b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_ASM; 13694b9ad928SBarry Smith pc->ops->setuponblocks = PCSetUpOnBlocks_ASM; 13704b9ad928SBarry Smith pc->ops->view = PCView_ASM; 13710a545947SLisandro Dalcin pc->ops->applyrichardson = NULL; 13724b9ad928SBarry Smith 1373bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetLocalSubdomains_C",PCASMSetLocalSubdomains_ASM);CHKERRQ(ierr); 1374bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetTotalSubdomains_C",PCASMSetTotalSubdomains_ASM);CHKERRQ(ierr); 1375bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetOverlap_C",PCASMSetOverlap_ASM);CHKERRQ(ierr); 1376bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetType_C",PCASMSetType_ASM);CHKERRQ(ierr); 1377c60c7ad4SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMGetType_C",PCASMGetType_ASM);CHKERRQ(ierr); 137812cd4985SMatthew G. Knepley ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetLocalType_C",PCASMSetLocalType_ASM);CHKERRQ(ierr); 137912cd4985SMatthew G. Knepley ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMGetLocalType_C",PCASMGetLocalType_ASM);CHKERRQ(ierr); 1380bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetSortIndices_C",PCASMSetSortIndices_ASM);CHKERRQ(ierr); 1381bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMGetSubKSP_C",PCASMGetSubKSP_ASM);CHKERRQ(ierr); 138280ec0b7dSPatrick Sanan ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMGetSubMatType_C",PCASMGetSubMatType_ASM);CHKERRQ(ierr); 138380ec0b7dSPatrick Sanan ierr = PetscObjectComposeFunction((PetscObject)pc,"PCASMSetSubMatType_C",PCASMSetSubMatType_ASM);CHKERRQ(ierr); 13844b9ad928SBarry Smith PetscFunctionReturn(0); 13854b9ad928SBarry Smith } 13864b9ad928SBarry Smith 138792bb6962SLisandro Dalcin /*@C 138892bb6962SLisandro Dalcin PCASMCreateSubdomains - Creates the index sets for the overlapping Schwarz 138992bb6962SLisandro Dalcin preconditioner for a any problem on a general grid. 139092bb6962SLisandro Dalcin 139192bb6962SLisandro Dalcin Collective 139292bb6962SLisandro Dalcin 139392bb6962SLisandro Dalcin Input Parameters: 139492bb6962SLisandro Dalcin + A - The global matrix operator 139592bb6962SLisandro Dalcin - n - the number of local blocks 139692bb6962SLisandro Dalcin 139792bb6962SLisandro Dalcin Output Parameters: 139892bb6962SLisandro Dalcin . outis - the array of index sets defining the subdomains 139992bb6962SLisandro Dalcin 140092bb6962SLisandro Dalcin Level: advanced 140192bb6962SLisandro Dalcin 14027d6bfa3bSBarry Smith Note: this generates nonoverlapping subdomains; the PCASM will generate the overlap 14037d6bfa3bSBarry Smith from these if you use PCASMSetLocalSubdomains() 14047d6bfa3bSBarry Smith 14057d6bfa3bSBarry Smith In the Fortran version you must provide the array outis[] already allocated of length n. 14067d6bfa3bSBarry Smith 140792bb6962SLisandro Dalcin .seealso: PCASMSetLocalSubdomains(), PCASMDestroySubdomains() 140892bb6962SLisandro Dalcin @*/ 14097087cfbeSBarry Smith PetscErrorCode PCASMCreateSubdomains(Mat A, PetscInt n, IS* outis[]) 141092bb6962SLisandro Dalcin { 141192bb6962SLisandro Dalcin MatPartitioning mpart; 141292bb6962SLisandro Dalcin const char *prefix; 141392bb6962SLisandro Dalcin PetscInt i,j,rstart,rend,bs; 1414976e8c5aSLisandro Dalcin PetscBool hasop, isbaij = PETSC_FALSE,foundpart = PETSC_FALSE; 14150298fd71SBarry Smith Mat Ad = NULL, adj; 141692bb6962SLisandro Dalcin IS ispart,isnumb,*is; 141792bb6962SLisandro Dalcin PetscErrorCode ierr; 141892bb6962SLisandro Dalcin 141992bb6962SLisandro Dalcin PetscFunctionBegin; 14200700a824SBarry Smith PetscValidHeaderSpecific(A,MAT_CLASSID,1); 142192bb6962SLisandro Dalcin PetscValidPointer(outis,3); 1422c1235816SBarry Smith if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of local blocks must be > 0, n = %D",n); 142392bb6962SLisandro Dalcin 142492bb6962SLisandro Dalcin /* Get prefix, row distribution, and block size */ 142592bb6962SLisandro Dalcin ierr = MatGetOptionsPrefix(A,&prefix);CHKERRQ(ierr); 142692bb6962SLisandro Dalcin ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 142792bb6962SLisandro Dalcin ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr); 142865e19b50SBarry 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); 142965e19b50SBarry Smith 143092bb6962SLisandro Dalcin /* Get diagonal block from matrix if possible */ 1431976e8c5aSLisandro Dalcin ierr = MatHasOperation(A,MATOP_GET_DIAGONAL_BLOCK,&hasop);CHKERRQ(ierr); 1432976e8c5aSLisandro Dalcin if (hasop) { 143311bd1e4dSLisandro Dalcin ierr = MatGetDiagonalBlock(A,&Ad);CHKERRQ(ierr); 143492bb6962SLisandro Dalcin } 143592bb6962SLisandro Dalcin if (Ad) { 1436b9e7e5c1SBarry Smith ierr = PetscObjectBaseTypeCompare((PetscObject)Ad,MATSEQBAIJ,&isbaij);CHKERRQ(ierr); 1437b9e7e5c1SBarry Smith if (!isbaij) {ierr = PetscObjectBaseTypeCompare((PetscObject)Ad,MATSEQSBAIJ,&isbaij);CHKERRQ(ierr);} 143892bb6962SLisandro Dalcin } 143992bb6962SLisandro Dalcin if (Ad && n > 1) { 1440ace3abfcSBarry Smith PetscBool match,done; 144192bb6962SLisandro Dalcin /* Try to setup a good matrix partitioning if available */ 144292bb6962SLisandro Dalcin ierr = MatPartitioningCreate(PETSC_COMM_SELF,&mpart);CHKERRQ(ierr); 144392bb6962SLisandro Dalcin ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr); 144492bb6962SLisandro Dalcin ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr); 1445251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGCURRENT,&match);CHKERRQ(ierr); 144692bb6962SLisandro Dalcin if (!match) { 1447251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGSQUARE,&match);CHKERRQ(ierr); 144892bb6962SLisandro Dalcin } 144992bb6962SLisandro Dalcin if (!match) { /* assume a "good" partitioner is available */ 14501a83f524SJed Brown PetscInt na; 14511a83f524SJed Brown const PetscInt *ia,*ja; 145292bb6962SLisandro Dalcin ierr = MatGetRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr); 145392bb6962SLisandro Dalcin if (done) { 145492bb6962SLisandro Dalcin /* Build adjacency matrix by hand. Unfortunately a call to 145592bb6962SLisandro Dalcin MatConvert(Ad,MATMPIADJ,MAT_INITIAL_MATRIX,&adj) will 145692bb6962SLisandro Dalcin remove the block-aij structure and we cannot expect 145792bb6962SLisandro Dalcin MatPartitioning to split vertices as we need */ 14580a545947SLisandro Dalcin PetscInt i,j,len,nnz,cnt,*iia=NULL,*jja=NULL; 14591a83f524SJed Brown const PetscInt *row; 146092bb6962SLisandro Dalcin nnz = 0; 146192bb6962SLisandro Dalcin for (i=0; i<na; i++) { /* count number of nonzeros */ 146292bb6962SLisandro Dalcin len = ia[i+1] - ia[i]; 146392bb6962SLisandro Dalcin row = ja + ia[i]; 146492bb6962SLisandro Dalcin for (j=0; j<len; j++) { 146592bb6962SLisandro Dalcin if (row[j] == i) { /* don't count diagonal */ 146692bb6962SLisandro Dalcin len--; break; 146792bb6962SLisandro Dalcin } 146892bb6962SLisandro Dalcin } 146992bb6962SLisandro Dalcin nnz += len; 147092bb6962SLisandro Dalcin } 1471854ce69bSBarry Smith ierr = PetscMalloc1(na+1,&iia);CHKERRQ(ierr); 1472854ce69bSBarry Smith ierr = PetscMalloc1(nnz,&jja);CHKERRQ(ierr); 147392bb6962SLisandro Dalcin nnz = 0; 147492bb6962SLisandro Dalcin iia[0] = 0; 147592bb6962SLisandro Dalcin for (i=0; i<na; i++) { /* fill adjacency */ 147692bb6962SLisandro Dalcin cnt = 0; 147792bb6962SLisandro Dalcin len = ia[i+1] - ia[i]; 147892bb6962SLisandro Dalcin row = ja + ia[i]; 147992bb6962SLisandro Dalcin for (j=0; j<len; j++) { 148092bb6962SLisandro Dalcin if (row[j] != i) { /* if not diagonal */ 148192bb6962SLisandro Dalcin jja[nnz+cnt++] = row[j]; 148292bb6962SLisandro Dalcin } 148392bb6962SLisandro Dalcin } 148492bb6962SLisandro Dalcin nnz += cnt; 148592bb6962SLisandro Dalcin iia[i+1] = nnz; 148692bb6962SLisandro Dalcin } 148792bb6962SLisandro Dalcin /* Partitioning of the adjacency matrix */ 14880298fd71SBarry Smith ierr = MatCreateMPIAdj(PETSC_COMM_SELF,na,na,iia,jja,NULL,&adj);CHKERRQ(ierr); 148992bb6962SLisandro Dalcin ierr = MatPartitioningSetAdjacency(mpart,adj);CHKERRQ(ierr); 149092bb6962SLisandro Dalcin ierr = MatPartitioningSetNParts(mpart,n);CHKERRQ(ierr); 149192bb6962SLisandro Dalcin ierr = MatPartitioningApply(mpart,&ispart);CHKERRQ(ierr); 149292bb6962SLisandro Dalcin ierr = ISPartitioningToNumbering(ispart,&isnumb);CHKERRQ(ierr); 1493fcfd50ebSBarry Smith ierr = MatDestroy(&adj);CHKERRQ(ierr); 149492bb6962SLisandro Dalcin foundpart = PETSC_TRUE; 149592bb6962SLisandro Dalcin } 149692bb6962SLisandro Dalcin ierr = MatRestoreRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr); 149792bb6962SLisandro Dalcin } 1498fcfd50ebSBarry Smith ierr = MatPartitioningDestroy(&mpart);CHKERRQ(ierr); 149992bb6962SLisandro Dalcin } 150092bb6962SLisandro Dalcin 1501785e854fSJed Brown ierr = PetscMalloc1(n,&is);CHKERRQ(ierr); 150292bb6962SLisandro Dalcin *outis = is; 150392bb6962SLisandro Dalcin 150492bb6962SLisandro Dalcin if (!foundpart) { 150592bb6962SLisandro Dalcin 150692bb6962SLisandro Dalcin /* Partitioning by contiguous chunks of rows */ 150792bb6962SLisandro Dalcin 150892bb6962SLisandro Dalcin PetscInt mbs = (rend-rstart)/bs; 150992bb6962SLisandro Dalcin PetscInt start = rstart; 151092bb6962SLisandro Dalcin for (i=0; i<n; i++) { 151192bb6962SLisandro Dalcin PetscInt count = (mbs/n + ((mbs % n) > i)) * bs; 151292bb6962SLisandro Dalcin ierr = ISCreateStride(PETSC_COMM_SELF,count,start,1,&is[i]);CHKERRQ(ierr); 151392bb6962SLisandro Dalcin start += count; 151492bb6962SLisandro Dalcin } 151592bb6962SLisandro Dalcin 151692bb6962SLisandro Dalcin } else { 151792bb6962SLisandro Dalcin 151892bb6962SLisandro Dalcin /* Partitioning by adjacency of diagonal block */ 151992bb6962SLisandro Dalcin 152092bb6962SLisandro Dalcin const PetscInt *numbering; 152192bb6962SLisandro Dalcin PetscInt *count,nidx,*indices,*newidx,start=0; 152292bb6962SLisandro Dalcin /* Get node count in each partition */ 1523785e854fSJed Brown ierr = PetscMalloc1(n,&count);CHKERRQ(ierr); 152492bb6962SLisandro Dalcin ierr = ISPartitioningCount(ispart,n,count);CHKERRQ(ierr); 152592bb6962SLisandro Dalcin if (isbaij && bs > 1) { /* adjust for the block-aij case */ 152692bb6962SLisandro Dalcin for (i=0; i<n; i++) count[i] *= bs; 152792bb6962SLisandro Dalcin } 152892bb6962SLisandro Dalcin /* Build indices from node numbering */ 152992bb6962SLisandro Dalcin ierr = ISGetLocalSize(isnumb,&nidx);CHKERRQ(ierr); 1530785e854fSJed Brown ierr = PetscMalloc1(nidx,&indices);CHKERRQ(ierr); 153192bb6962SLisandro Dalcin for (i=0; i<nidx; i++) indices[i] = i; /* needs to be initialized */ 153292bb6962SLisandro Dalcin ierr = ISGetIndices(isnumb,&numbering);CHKERRQ(ierr); 153392bb6962SLisandro Dalcin ierr = PetscSortIntWithPermutation(nidx,numbering,indices);CHKERRQ(ierr); 153492bb6962SLisandro Dalcin ierr = ISRestoreIndices(isnumb,&numbering);CHKERRQ(ierr); 153592bb6962SLisandro Dalcin if (isbaij && bs > 1) { /* adjust for the block-aij case */ 1536785e854fSJed Brown ierr = PetscMalloc1(nidx*bs,&newidx);CHKERRQ(ierr); 15372fa5cd67SKarl Rupp for (i=0; i<nidx; i++) { 15382fa5cd67SKarl Rupp for (j=0; j<bs; j++) newidx[i*bs+j] = indices[i]*bs + j; 15392fa5cd67SKarl Rupp } 154092bb6962SLisandro Dalcin ierr = PetscFree(indices);CHKERRQ(ierr); 154192bb6962SLisandro Dalcin nidx *= bs; 154292bb6962SLisandro Dalcin indices = newidx; 154392bb6962SLisandro Dalcin } 154492bb6962SLisandro Dalcin /* Shift to get global indices */ 154592bb6962SLisandro Dalcin for (i=0; i<nidx; i++) indices[i] += rstart; 154692bb6962SLisandro Dalcin 154792bb6962SLisandro Dalcin /* Build the index sets for each block */ 154892bb6962SLisandro Dalcin for (i=0; i<n; i++) { 154970b3c8c7SBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,count[i],&indices[start],PETSC_COPY_VALUES,&is[i]);CHKERRQ(ierr); 155092bb6962SLisandro Dalcin ierr = ISSort(is[i]);CHKERRQ(ierr); 155192bb6962SLisandro Dalcin start += count[i]; 155292bb6962SLisandro Dalcin } 155392bb6962SLisandro Dalcin 15543bf036e2SBarry Smith ierr = PetscFree(count);CHKERRQ(ierr); 15553bf036e2SBarry Smith ierr = PetscFree(indices);CHKERRQ(ierr); 1556fcfd50ebSBarry Smith ierr = ISDestroy(&isnumb);CHKERRQ(ierr); 1557fcfd50ebSBarry Smith ierr = ISDestroy(&ispart);CHKERRQ(ierr); 155892bb6962SLisandro Dalcin 155992bb6962SLisandro Dalcin } 156092bb6962SLisandro Dalcin PetscFunctionReturn(0); 156192bb6962SLisandro Dalcin } 156292bb6962SLisandro Dalcin 156392bb6962SLisandro Dalcin /*@C 156492bb6962SLisandro Dalcin PCASMDestroySubdomains - Destroys the index sets created with 156592bb6962SLisandro Dalcin PCASMCreateSubdomains(). Should be called after setting subdomains 156692bb6962SLisandro Dalcin with PCASMSetLocalSubdomains(). 156792bb6962SLisandro Dalcin 156892bb6962SLisandro Dalcin Collective 156992bb6962SLisandro Dalcin 157092bb6962SLisandro Dalcin Input Parameters: 157192bb6962SLisandro Dalcin + n - the number of index sets 15722b691e39SMatthew Knepley . is - the array of index sets 15730298fd71SBarry Smith - is_local - the array of local index sets, can be NULL 157492bb6962SLisandro Dalcin 157592bb6962SLisandro Dalcin Level: advanced 157692bb6962SLisandro Dalcin 157792bb6962SLisandro Dalcin .seealso: PCASMCreateSubdomains(), PCASMSetLocalSubdomains() 157892bb6962SLisandro Dalcin @*/ 15797087cfbeSBarry Smith PetscErrorCode PCASMDestroySubdomains(PetscInt n, IS is[], IS is_local[]) 158092bb6962SLisandro Dalcin { 158192bb6962SLisandro Dalcin PetscInt i; 158292bb6962SLisandro Dalcin PetscErrorCode ierr; 15835fd66863SKarl Rupp 158492bb6962SLisandro Dalcin PetscFunctionBegin; 1585a35b7b57SDmitry Karpeev if (n <= 0) PetscFunctionReturn(0); 1586a35b7b57SDmitry Karpeev if (is) { 158792bb6962SLisandro Dalcin PetscValidPointer(is,2); 1588fcfd50ebSBarry Smith for (i=0; i<n; i++) { ierr = ISDestroy(&is[i]);CHKERRQ(ierr); } 158992bb6962SLisandro Dalcin ierr = PetscFree(is);CHKERRQ(ierr); 1590a35b7b57SDmitry Karpeev } 15912b691e39SMatthew Knepley if (is_local) { 15922b691e39SMatthew Knepley PetscValidPointer(is_local,3); 1593fcfd50ebSBarry Smith for (i=0; i<n; i++) { ierr = ISDestroy(&is_local[i]);CHKERRQ(ierr); } 15942b691e39SMatthew Knepley ierr = PetscFree(is_local);CHKERRQ(ierr); 15952b691e39SMatthew Knepley } 159692bb6962SLisandro Dalcin PetscFunctionReturn(0); 159792bb6962SLisandro Dalcin } 159892bb6962SLisandro Dalcin 15994b9ad928SBarry Smith /*@ 16004b9ad928SBarry Smith PCASMCreateSubdomains2D - Creates the index sets for the overlapping Schwarz 16014b9ad928SBarry Smith preconditioner for a two-dimensional problem on a regular grid. 16024b9ad928SBarry Smith 16034b9ad928SBarry Smith Not Collective 16044b9ad928SBarry Smith 16054b9ad928SBarry Smith Input Parameters: 16064b9ad928SBarry Smith + m, n - the number of mesh points in the x and y directions 16074b9ad928SBarry Smith . M, N - the number of subdomains in the x and y directions 16084b9ad928SBarry Smith . dof - degrees of freedom per node 16094b9ad928SBarry Smith - overlap - overlap in mesh lines 16104b9ad928SBarry Smith 16114b9ad928SBarry Smith Output Parameters: 16124b9ad928SBarry Smith + Nsub - the number of subdomains created 16133d03bb48SJed Brown . is - array of index sets defining overlapping (if overlap > 0) subdomains 16143d03bb48SJed Brown - is_local - array of index sets defining non-overlapping subdomains 16154b9ad928SBarry Smith 16164b9ad928SBarry Smith Note: 16174b9ad928SBarry Smith Presently PCAMSCreateSubdomains2d() is valid only for sequential 16184b9ad928SBarry Smith preconditioners. More general related routines are 16194b9ad928SBarry Smith PCASMSetTotalSubdomains() and PCASMSetLocalSubdomains(). 16204b9ad928SBarry Smith 16214b9ad928SBarry Smith Level: advanced 16224b9ad928SBarry Smith 16234b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetLocalSubdomains(), PCASMGetSubKSP(), 16244b9ad928SBarry Smith PCASMSetOverlap() 16254b9ad928SBarry Smith @*/ 16267087cfbeSBarry Smith PetscErrorCode PCASMCreateSubdomains2D(PetscInt m,PetscInt n,PetscInt M,PetscInt N,PetscInt dof,PetscInt overlap,PetscInt *Nsub,IS **is,IS **is_local) 16274b9ad928SBarry Smith { 16283d03bb48SJed Brown PetscInt i,j,height,width,ystart,xstart,yleft,yright,xleft,xright,loc_outer; 16296849ba73SBarry Smith PetscErrorCode ierr; 163013f74950SBarry Smith PetscInt nidx,*idx,loc,ii,jj,count; 16314b9ad928SBarry Smith 16324b9ad928SBarry Smith PetscFunctionBegin; 1633e32f2f54SBarry Smith if (dof != 1) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_SUP," "); 16344b9ad928SBarry Smith 16354b9ad928SBarry Smith *Nsub = N*M; 1636854ce69bSBarry Smith ierr = PetscMalloc1(*Nsub,is);CHKERRQ(ierr); 1637854ce69bSBarry Smith ierr = PetscMalloc1(*Nsub,is_local);CHKERRQ(ierr); 16384b9ad928SBarry Smith ystart = 0; 16393d03bb48SJed Brown loc_outer = 0; 16404b9ad928SBarry Smith for (i=0; i<N; i++) { 16414b9ad928SBarry Smith height = n/N + ((n % N) > i); /* height of subdomain */ 1642e32f2f54SBarry Smith if (height < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many N subdomains for mesh dimension n"); 16434b9ad928SBarry Smith yleft = ystart - overlap; if (yleft < 0) yleft = 0; 16444b9ad928SBarry Smith yright = ystart + height + overlap; if (yright > n) yright = n; 16454b9ad928SBarry Smith xstart = 0; 16464b9ad928SBarry Smith for (j=0; j<M; j++) { 16474b9ad928SBarry Smith width = m/M + ((m % M) > j); /* width of subdomain */ 1648e32f2f54SBarry Smith if (width < 2) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many M subdomains for mesh dimension m"); 16494b9ad928SBarry Smith xleft = xstart - overlap; if (xleft < 0) xleft = 0; 16504b9ad928SBarry Smith xright = xstart + width + overlap; if (xright > m) xright = m; 16514b9ad928SBarry Smith nidx = (xright - xleft)*(yright - yleft); 1652785e854fSJed Brown ierr = PetscMalloc1(nidx,&idx);CHKERRQ(ierr); 16534b9ad928SBarry Smith loc = 0; 16544b9ad928SBarry Smith for (ii=yleft; ii<yright; ii++) { 16554b9ad928SBarry Smith count = m*ii + xleft; 16562fa5cd67SKarl Rupp for (jj=xleft; jj<xright; jj++) idx[loc++] = count++; 16574b9ad928SBarry Smith } 165870b3c8c7SBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,nidx,idx,PETSC_COPY_VALUES,(*is)+loc_outer);CHKERRQ(ierr); 16593d03bb48SJed Brown if (overlap == 0) { 16603d03bb48SJed Brown ierr = PetscObjectReference((PetscObject)(*is)[loc_outer]);CHKERRQ(ierr); 16612fa5cd67SKarl Rupp 16623d03bb48SJed Brown (*is_local)[loc_outer] = (*is)[loc_outer]; 16633d03bb48SJed Brown } else { 16643d03bb48SJed Brown for (loc=0,ii=ystart; ii<ystart+height; ii++) { 16653d03bb48SJed Brown for (jj=xstart; jj<xstart+width; jj++) { 16663d03bb48SJed Brown idx[loc++] = m*ii + jj; 16673d03bb48SJed Brown } 16683d03bb48SJed Brown } 166970b3c8c7SBarry Smith ierr = ISCreateGeneral(PETSC_COMM_SELF,loc,idx,PETSC_COPY_VALUES,*is_local+loc_outer);CHKERRQ(ierr); 16703d03bb48SJed Brown } 16714b9ad928SBarry Smith ierr = PetscFree(idx);CHKERRQ(ierr); 16724b9ad928SBarry Smith xstart += width; 16733d03bb48SJed Brown loc_outer++; 16744b9ad928SBarry Smith } 16754b9ad928SBarry Smith ystart += height; 16764b9ad928SBarry Smith } 16774b9ad928SBarry Smith for (i=0; i<*Nsub; i++) { ierr = ISSort((*is)[i]);CHKERRQ(ierr); } 16784b9ad928SBarry Smith PetscFunctionReturn(0); 16794b9ad928SBarry Smith } 16804b9ad928SBarry Smith 16814b9ad928SBarry Smith /*@C 16824b9ad928SBarry Smith PCASMGetLocalSubdomains - Gets the local subdomains (for this processor 16834b9ad928SBarry Smith only) for the additive Schwarz preconditioner. 16844b9ad928SBarry Smith 1685ad4df100SBarry Smith Not Collective 16864b9ad928SBarry Smith 16874b9ad928SBarry Smith Input Parameter: 16884b9ad928SBarry Smith . pc - the preconditioner context 16894b9ad928SBarry Smith 16904b9ad928SBarry Smith Output Parameters: 16914b9ad928SBarry Smith + n - the number of subdomains for this processor (default value = 1) 16922b691e39SMatthew Knepley . is - the index sets that define the subdomains for this processor 16930298fd71SBarry Smith - is_local - the index sets that define the local part of the subdomains for this processor (can be NULL) 16944b9ad928SBarry Smith 16954b9ad928SBarry Smith 16964b9ad928SBarry Smith Notes: 16974b9ad928SBarry Smith The IS numbering is in the parallel, global numbering of the vector. 16984b9ad928SBarry Smith 16994b9ad928SBarry Smith Level: advanced 17004b9ad928SBarry Smith 17014b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(), 17024b9ad928SBarry Smith PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubmatrices() 17034b9ad928SBarry Smith @*/ 17047087cfbeSBarry Smith PetscErrorCode PCASMGetLocalSubdomains(PC pc,PetscInt *n,IS *is[],IS *is_local[]) 17054b9ad928SBarry Smith { 17062a808120SBarry Smith PC_ASM *osm = (PC_ASM*)pc->data; 170792bb6962SLisandro Dalcin PetscErrorCode ierr; 1708ace3abfcSBarry Smith PetscBool match; 17094b9ad928SBarry Smith 17104b9ad928SBarry Smith PetscFunctionBegin; 17110700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 17124482741eSBarry Smith PetscValidIntPointer(n,2); 171392bb6962SLisandro Dalcin if (is) PetscValidPointer(is,3); 1714251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)pc,PCASM,&match);CHKERRQ(ierr); 17152a808120SBarry Smith if (!match) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG,"PC is not a PCASM"); 17164b9ad928SBarry Smith if (n) *n = osm->n_local_true; 17174b9ad928SBarry Smith if (is) *is = osm->is; 17182b691e39SMatthew Knepley if (is_local) *is_local = osm->is_local; 17194b9ad928SBarry Smith PetscFunctionReturn(0); 17204b9ad928SBarry Smith } 17214b9ad928SBarry Smith 17224b9ad928SBarry Smith /*@C 17234b9ad928SBarry Smith PCASMGetLocalSubmatrices - Gets the local submatrices (for this processor 17244b9ad928SBarry Smith only) for the additive Schwarz preconditioner. 17254b9ad928SBarry Smith 1726ad4df100SBarry Smith Not Collective 17274b9ad928SBarry Smith 17284b9ad928SBarry Smith Input Parameter: 17294b9ad928SBarry Smith . pc - the preconditioner context 17304b9ad928SBarry Smith 17314b9ad928SBarry Smith Output Parameters: 17324b9ad928SBarry Smith + n - the number of matrices for this processor (default value = 1) 17334b9ad928SBarry Smith - mat - the matrices 17344b9ad928SBarry Smith 17354b9ad928SBarry Smith Level: advanced 17364b9ad928SBarry Smith 173795452b02SPatrick Sanan Notes: 173834a84908Sprj- Call after PCSetUp() (or KSPSetUp()) but before PCApply() and before PCSetUpOnBlocks()) 1739cf739d55SBarry Smith 174034a84908Sprj- Usually one would use PCSetModifySubMatrices() to change the submatrices in building the preconditioner. 1741cf739d55SBarry Smith 17424b9ad928SBarry Smith .seealso: PCASMSetTotalSubdomains(), PCASMSetOverlap(), PCASMGetSubKSP(), 174334a84908Sprj- PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubdomains(), PCSetModifySubMatrices() 17444b9ad928SBarry Smith @*/ 17457087cfbeSBarry Smith PetscErrorCode PCASMGetLocalSubmatrices(PC pc,PetscInt *n,Mat *mat[]) 17464b9ad928SBarry Smith { 17474b9ad928SBarry Smith PC_ASM *osm; 174892bb6962SLisandro Dalcin PetscErrorCode ierr; 1749ace3abfcSBarry Smith PetscBool match; 17504b9ad928SBarry Smith 17514b9ad928SBarry Smith PetscFunctionBegin; 17520700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 175392bb6962SLisandro Dalcin PetscValidIntPointer(n,2); 175492bb6962SLisandro Dalcin if (mat) PetscValidPointer(mat,3); 175534a84908Sprj- if (!pc->setupcalled) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call after KSPSetUp() or PCSetUp()."); 1756251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)pc,PCASM,&match);CHKERRQ(ierr); 175792bb6962SLisandro Dalcin if (!match) { 175892bb6962SLisandro Dalcin if (n) *n = 0; 17590298fd71SBarry Smith if (mat) *mat = NULL; 176092bb6962SLisandro Dalcin } else { 17614b9ad928SBarry Smith osm = (PC_ASM*)pc->data; 17624b9ad928SBarry Smith if (n) *n = osm->n_local_true; 17634b9ad928SBarry Smith if (mat) *mat = osm->pmat; 176492bb6962SLisandro Dalcin } 17654b9ad928SBarry Smith PetscFunctionReturn(0); 17664b9ad928SBarry Smith } 1767d709ea83SDmitry Karpeev 1768d709ea83SDmitry Karpeev /*@ 1769d709ea83SDmitry Karpeev PCASMSetDMSubdomains - Indicates whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible. 1770f1ee410cSBarry Smith 1771d709ea83SDmitry Karpeev Logically Collective 1772d709ea83SDmitry Karpeev 1773d709ea83SDmitry Karpeev Input Parameter: 1774d709ea83SDmitry Karpeev + pc - the preconditioner 1775d709ea83SDmitry Karpeev - flg - boolean indicating whether to use subdomains defined by the DM 1776d709ea83SDmitry Karpeev 1777d709ea83SDmitry Karpeev Options Database Key: 1778d709ea83SDmitry Karpeev . -pc_asm_dm_subdomains 1779d709ea83SDmitry Karpeev 1780d709ea83SDmitry Karpeev Level: intermediate 1781d709ea83SDmitry Karpeev 1782d709ea83SDmitry Karpeev Notes: 1783d709ea83SDmitry Karpeev PCASMSetTotalSubdomains() and PCASMSetOverlap() take precedence over PCASMSetDMSubdomains(), 1784d709ea83SDmitry Karpeev so setting either of the first two effectively turns the latter off. 1785d709ea83SDmitry Karpeev 1786d709ea83SDmitry Karpeev .seealso: PCASMGetDMSubdomains(), PCASMSetTotalSubdomains(), PCASMSetOverlap() 1787d709ea83SDmitry Karpeev PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubdomains() 1788d709ea83SDmitry Karpeev @*/ 1789d709ea83SDmitry Karpeev PetscErrorCode PCASMSetDMSubdomains(PC pc,PetscBool flg) 1790d709ea83SDmitry Karpeev { 1791d709ea83SDmitry Karpeev PC_ASM *osm = (PC_ASM*)pc->data; 1792d709ea83SDmitry Karpeev PetscErrorCode ierr; 1793d709ea83SDmitry Karpeev PetscBool match; 1794d709ea83SDmitry Karpeev 1795d709ea83SDmitry Karpeev PetscFunctionBegin; 1796d709ea83SDmitry Karpeev PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1797d709ea83SDmitry Karpeev PetscValidLogicalCollectiveBool(pc,flg,2); 1798d709ea83SDmitry Karpeev if (pc->setupcalled) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Not for a setup PC."); 1799d709ea83SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)pc,PCASM,&match);CHKERRQ(ierr); 1800d709ea83SDmitry Karpeev if (match) { 1801d709ea83SDmitry Karpeev osm->dm_subdomains = flg; 1802d709ea83SDmitry Karpeev } 1803d709ea83SDmitry Karpeev PetscFunctionReturn(0); 1804d709ea83SDmitry Karpeev } 1805d709ea83SDmitry Karpeev 1806d709ea83SDmitry Karpeev /*@ 1807d709ea83SDmitry Karpeev PCASMGetDMSubdomains - Returns flag indicating whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible. 1808d709ea83SDmitry Karpeev Not Collective 1809d709ea83SDmitry Karpeev 1810d709ea83SDmitry Karpeev Input Parameter: 1811d709ea83SDmitry Karpeev . pc - the preconditioner 1812d709ea83SDmitry Karpeev 1813d709ea83SDmitry Karpeev Output Parameter: 1814d709ea83SDmitry Karpeev . flg - boolean indicating whether to use subdomains defined by the DM 1815d709ea83SDmitry Karpeev 1816d709ea83SDmitry Karpeev Level: intermediate 1817d709ea83SDmitry Karpeev 1818d709ea83SDmitry Karpeev .seealso: PCASMSetDMSubdomains(), PCASMSetTotalSubdomains(), PCASMSetOverlap() 1819d709ea83SDmitry Karpeev PCASMCreateSubdomains2D(), PCASMSetLocalSubdomains(), PCASMGetLocalSubdomains() 1820d709ea83SDmitry Karpeev @*/ 1821d709ea83SDmitry Karpeev PetscErrorCode PCASMGetDMSubdomains(PC pc,PetscBool* flg) 1822d709ea83SDmitry Karpeev { 1823d709ea83SDmitry Karpeev PC_ASM *osm = (PC_ASM*)pc->data; 1824d709ea83SDmitry Karpeev PetscErrorCode ierr; 1825d709ea83SDmitry Karpeev PetscBool match; 1826d709ea83SDmitry Karpeev 1827d709ea83SDmitry Karpeev PetscFunctionBegin; 1828d709ea83SDmitry Karpeev PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1829534a8f05SLisandro Dalcin PetscValidBoolPointer(flg,2); 1830d709ea83SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)pc,PCASM,&match);CHKERRQ(ierr); 183156ea09ceSStefano Zampini if (match) *flg = osm->dm_subdomains; 183256ea09ceSStefano Zampini else *flg = PETSC_FALSE; 1833d709ea83SDmitry Karpeev PetscFunctionReturn(0); 1834d709ea83SDmitry Karpeev } 183580ec0b7dSPatrick Sanan 183680ec0b7dSPatrick Sanan /*@ 183780ec0b7dSPatrick Sanan PCASMGetSubMatType - Gets the matrix type used for ASM subsolves, as a string. 183880ec0b7dSPatrick Sanan 183980ec0b7dSPatrick Sanan Not Collective 184080ec0b7dSPatrick Sanan 184180ec0b7dSPatrick Sanan Input Parameter: 184280ec0b7dSPatrick Sanan . pc - the PC 184380ec0b7dSPatrick Sanan 184480ec0b7dSPatrick Sanan Output Parameter: 1845f1ee410cSBarry Smith . -pc_asm_sub_mat_type - name of matrix type 184680ec0b7dSPatrick Sanan 184780ec0b7dSPatrick Sanan Level: advanced 184880ec0b7dSPatrick Sanan 184980ec0b7dSPatrick Sanan .seealso: PCASMSetSubMatType(), PCASM, PCSetType(), VecSetType(), MatType, Mat 185080ec0b7dSPatrick Sanan @*/ 185156ea09ceSStefano Zampini PetscErrorCode PCASMGetSubMatType(PC pc,MatType *sub_mat_type) 185256ea09ceSStefano Zampini { 185380ec0b7dSPatrick Sanan PetscErrorCode ierr; 185480ec0b7dSPatrick Sanan 185556ea09ceSStefano Zampini PetscFunctionBegin; 185656ea09ceSStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 185780ec0b7dSPatrick Sanan ierr = PetscTryMethod(pc,"PCASMGetSubMatType_C",(PC,MatType*),(pc,sub_mat_type));CHKERRQ(ierr); 185880ec0b7dSPatrick Sanan PetscFunctionReturn(0); 185980ec0b7dSPatrick Sanan } 186080ec0b7dSPatrick Sanan 186180ec0b7dSPatrick Sanan /*@ 186280ec0b7dSPatrick Sanan PCASMSetSubMatType - Set the type of matrix used for ASM subsolves 186380ec0b7dSPatrick Sanan 186480ec0b7dSPatrick Sanan Collective on Mat 186580ec0b7dSPatrick Sanan 186680ec0b7dSPatrick Sanan Input Parameters: 186780ec0b7dSPatrick Sanan + pc - the PC object 186880ec0b7dSPatrick Sanan - sub_mat_type - matrix type 186980ec0b7dSPatrick Sanan 187080ec0b7dSPatrick Sanan Options Database Key: 187180ec0b7dSPatrick 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. 187280ec0b7dSPatrick Sanan 187380ec0b7dSPatrick Sanan Notes: 187480ec0b7dSPatrick Sanan See "${PETSC_DIR}/include/petscmat.h" for available types 187580ec0b7dSPatrick Sanan 187680ec0b7dSPatrick Sanan Level: advanced 187780ec0b7dSPatrick Sanan 187880ec0b7dSPatrick Sanan .seealso: PCASMGetSubMatType(), PCASM, PCSetType(), VecSetType(), MatType, Mat 187980ec0b7dSPatrick Sanan @*/ 188080ec0b7dSPatrick Sanan PetscErrorCode PCASMSetSubMatType(PC pc,MatType sub_mat_type) 188180ec0b7dSPatrick Sanan { 188280ec0b7dSPatrick Sanan PetscErrorCode ierr; 188380ec0b7dSPatrick Sanan 188456ea09ceSStefano Zampini PetscFunctionBegin; 188556ea09ceSStefano Zampini PetscValidHeaderSpecific(pc,PC_CLASSID,1); 188680ec0b7dSPatrick Sanan ierr = PetscTryMethod(pc,"PCASMSetSubMatType_C",(PC,MatType),(pc,sub_mat_type));CHKERRQ(ierr); 188780ec0b7dSPatrick Sanan PetscFunctionReturn(0); 188880ec0b7dSPatrick Sanan } 1889