1b1a0a8a3SJed Brown /* 2f746d493SDmitry Karpeev This file defines an "generalized" additive Schwarz preconditioner for any Mat implementation. 36a4f0f73SDmitry Karpeev In this version each processor may intersect multiple subdomains and any subdomain may 46a4f0f73SDmitry Karpeev intersect multiple processors. Intersections of subdomains with processors are called *local 56a4f0f73SDmitry Karpeev subdomains*. 6b1a0a8a3SJed Brown 78f3b4b4dSDmitry Karpeev N - total number of distinct global subdomains (set explicitly in PCGASMSetTotalSubdomains() or implicitly PCGASMSetSubdomains() and then calculated in PCSetUp_GASM()) 86a4f0f73SDmitry Karpeev n - actual number of local subdomains on this processor (set in PCGASMSetSubdomains() or calculated in PCGASMSetTotalSubdomains()) 98f3b4b4dSDmitry Karpeev nmax - maximum number of local subdomains per processor (calculated in PCSetUp_GASM()) 10b1a0a8a3SJed Brown */ 11af0996ceSBarry Smith #include <petsc/private/pcimpl.h> /*I "petscpc.h" I*/ 121e25c274SJed Brown #include <petscdm.h> 13b1a0a8a3SJed Brown 14b1a0a8a3SJed Brown typedef struct { 15f746d493SDmitry Karpeev PetscInt N,n,nmax; 16b1a0a8a3SJed Brown PetscInt overlap; /* overlap requested by user */ 178f3b4b4dSDmitry Karpeev PCGASMType type; /* use reduced interpolation, restriction or both */ 188f3b4b4dSDmitry Karpeev PetscBool type_set; /* if user set this value (so won't change it for symmetric problems) */ 198f3b4b4dSDmitry Karpeev PetscBool same_subdomain_solvers; /* flag indicating whether all local solvers are same */ 208f3b4b4dSDmitry Karpeev PetscBool sort_indices; /* flag to sort subdomain indices */ 218f3b4b4dSDmitry Karpeev PetscBool user_subdomains; /* whether the user set explicit subdomain index sets -- keep them on PCReset() */ 228f3b4b4dSDmitry Karpeev PetscBool dm_subdomains; /* whether DM is allowed to define subdomains */ 238f3b4b4dSDmitry Karpeev IS *ois; /* index sets that define the outer (conceptually, overlapping) subdomains */ 248f3b4b4dSDmitry Karpeev IS *iis; /* index sets that define the inner (conceptually, nonoverlapping) subdomains */ 258f3b4b4dSDmitry Karpeev KSP *ksp; /* linear solvers for each subdomain */ 268f3b4b4dSDmitry Karpeev Mat *pmat; /* subdomain block matrices */ 27f746d493SDmitry Karpeev Vec gx,gy; /* Merged work vectors */ 28f746d493SDmitry Karpeev Vec *x,*y; /* Split work vectors; storage aliases pieces of storage of the above merged vectors. */ 296a4f0f73SDmitry Karpeev VecScatter gorestriction; /* merged restriction to disjoint union of outer subdomains */ 306a4f0f73SDmitry Karpeev VecScatter girestriction; /* merged restriction to disjoint union of inner subdomains */ 31f746d493SDmitry Karpeev } PC_GASM; 32b1a0a8a3SJed Brown 33b1a0a8a3SJed Brown #undef __FUNCT__ 348f3b4b4dSDmitry Karpeev #define __FUNCT__ "PCGASMComputeGlobalSubdomainNumbering_Private" 358f3b4b4dSDmitry Karpeev static PetscErrorCode PCGASMComputeGlobalSubdomainNumbering_Private(PC pc,PetscInt **numbering,PetscInt **permutation) 368f3b4b4dSDmitry Karpeev { 378f3b4b4dSDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 388f3b4b4dSDmitry Karpeev PetscInt i; 398f3b4b4dSDmitry Karpeev PetscErrorCode ierr; 408f3b4b4dSDmitry Karpeev 418f3b4b4dSDmitry Karpeev PetscFunctionBegin; 428f3b4b4dSDmitry Karpeev /* Determine the number of globally-distinct subdomains and compute a global numbering for them. */ 438f3b4b4dSDmitry Karpeev ierr = PetscMalloc2(osm->n,numbering,osm->n,permutation);CHKERRQ(ierr); 448f3b4b4dSDmitry Karpeev ierr = PetscObjectsListGetGlobalNumbering(PetscObjectComm((PetscObject)pc),osm->n,(PetscObject*)osm->iis,NULL,*numbering);CHKERRQ(ierr); 458f3b4b4dSDmitry Karpeev for (i = 0; i < osm->n; ++i) (*permutation)[i] = i; 468f3b4b4dSDmitry Karpeev ierr = PetscSortIntWithPermutation(osm->n,*numbering,*permutation);CHKERRQ(ierr); 478f3b4b4dSDmitry Karpeev PetscFunctionReturn(0); 488f3b4b4dSDmitry Karpeev } 498f3b4b4dSDmitry Karpeev 508f3b4b4dSDmitry Karpeev #undef __FUNCT__ 5106b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMSubdomainView_Private" 5206b43e7eSDmitry Karpeev static PetscErrorCode PCGASMSubdomainView_Private(PC pc, PetscInt i, PetscViewer viewer) 53af538404SDmitry Karpeev { 54af538404SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 5506b43e7eSDmitry Karpeev PetscInt j,nidx; 56af538404SDmitry Karpeev const PetscInt *idx; 5706b43e7eSDmitry Karpeev PetscViewer sviewer; 58af538404SDmitry Karpeev char *cidx; 59af538404SDmitry Karpeev PetscErrorCode ierr; 60af538404SDmitry Karpeev 61af538404SDmitry Karpeev PetscFunctionBegin; 62ce94432eSBarry Smith if (i < 0 || i > osm->n) SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONG, "Invalid subdomain %D: must nonnegative and less than %D", i, osm->n); 634bde246dSDmitry Karpeev /* Inner subdomains. */ 644bde246dSDmitry Karpeev ierr = ISGetLocalSize(osm->iis[i], &nidx);CHKERRQ(ierr); 654bde246dSDmitry Karpeev /* 664bde246dSDmitry Karpeev No more than 15 characters per index plus a space. 674bde246dSDmitry Karpeev PetscViewerStringSPrintf requires a string of size at least 2, so use (nidx+1) instead of nidx, 684bde246dSDmitry Karpeev in case nidx == 0. That will take care of the space for the trailing '\0' as well. 694bde246dSDmitry Karpeev For nidx == 0, the whole string 16 '\0'. 704bde246dSDmitry Karpeev */ 7189d949e2SBarry Smith ierr = PetscMalloc1(16*(nidx+1)+1, &cidx);CHKERRQ(ierr); 724bde246dSDmitry Karpeev ierr = PetscViewerStringOpen(PETSC_COMM_SELF, cidx, 16*(nidx+1)+1, &sviewer);CHKERRQ(ierr); 734bde246dSDmitry Karpeev ierr = ISGetIndices(osm->iis[i], &idx);CHKERRQ(ierr); 744bde246dSDmitry Karpeev for (j = 0; j < nidx; ++j) { 754bde246dSDmitry Karpeev ierr = PetscViewerStringSPrintf(sviewer, "%D ", idx[j]);CHKERRQ(ierr); 764bde246dSDmitry Karpeev } 774bde246dSDmitry Karpeev ierr = ISRestoreIndices(osm->iis[i],&idx);CHKERRQ(ierr); 784bde246dSDmitry Karpeev ierr = PetscViewerDestroy(&sviewer);CHKERRQ(ierr); 794bde246dSDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer, "Inner subdomain:\n");CHKERRQ(ierr); 804bde246dSDmitry Karpeev ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 811575c14dSBarry Smith ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 824bde246dSDmitry Karpeev ierr = PetscViewerASCIISynchronizedPrintf(viewer, "%s", cidx);CHKERRQ(ierr); 834bde246dSDmitry Karpeev ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 841575c14dSBarry Smith ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 854bde246dSDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer, "\n");CHKERRQ(ierr); 864bde246dSDmitry Karpeev ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 874bde246dSDmitry Karpeev ierr = PetscFree(cidx);CHKERRQ(ierr); 884bde246dSDmitry Karpeev /* Outer subdomains. */ 896a4f0f73SDmitry Karpeev ierr = ISGetLocalSize(osm->ois[i], &nidx);CHKERRQ(ierr); 90eec7e3faSDmitry Karpeev /* 91eec7e3faSDmitry Karpeev No more than 15 characters per index plus a space. 92eec7e3faSDmitry Karpeev PetscViewerStringSPrintf requires a string of size at least 2, so use (nidx+1) instead of nidx, 93eec7e3faSDmitry Karpeev in case nidx == 0. That will take care of the space for the trailing '\0' as well. 94eec7e3faSDmitry Karpeev For nidx == 0, the whole string 16 '\0'. 95eec7e3faSDmitry Karpeev */ 9689d949e2SBarry Smith ierr = PetscMalloc1(16*(nidx+1)+1, &cidx);CHKERRQ(ierr); 9706b43e7eSDmitry Karpeev ierr = PetscViewerStringOpen(PETSC_COMM_SELF, cidx, 16*(nidx+1)+1, &sviewer);CHKERRQ(ierr); 986a4f0f73SDmitry Karpeev ierr = ISGetIndices(osm->ois[i], &idx);CHKERRQ(ierr); 99af538404SDmitry Karpeev for (j = 0; j < nidx; ++j) { 100af538404SDmitry Karpeev ierr = PetscViewerStringSPrintf(sviewer,"%D ", idx[j]);CHKERRQ(ierr); 101af538404SDmitry Karpeev } 1026bf464f9SBarry Smith ierr = PetscViewerDestroy(&sviewer);CHKERRQ(ierr); 1036a4f0f73SDmitry Karpeev ierr = ISRestoreIndices(osm->ois[i],&idx);CHKERRQ(ierr); 1046a4f0f73SDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer, "Outer subdomain:\n");CHKERRQ(ierr); 10506b43e7eSDmitry Karpeev ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 1061575c14dSBarry Smith ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 107af538404SDmitry Karpeev ierr = PetscViewerASCIISynchronizedPrintf(viewer, "%s", cidx);CHKERRQ(ierr); 108af538404SDmitry Karpeev ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 1091575c14dSBarry Smith ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 110af538404SDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer, "\n");CHKERRQ(ierr); 11106b43e7eSDmitry Karpeev ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 112af538404SDmitry Karpeev ierr = PetscFree(cidx);CHKERRQ(ierr); 11306b43e7eSDmitry Karpeev PetscFunctionReturn(0); 11406b43e7eSDmitry Karpeev } 11506b43e7eSDmitry Karpeev 11606b43e7eSDmitry Karpeev #undef __FUNCT__ 11706b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMPrintSubdomains" 11806b43e7eSDmitry Karpeev static PetscErrorCode PCGASMPrintSubdomains(PC pc) 11906b43e7eSDmitry Karpeev { 12006b43e7eSDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 12106b43e7eSDmitry Karpeev const char *prefix; 12206b43e7eSDmitry Karpeev char fname[PETSC_MAX_PATH_LEN+1]; 1238f3b4b4dSDmitry Karpeev PetscInt l, d, count; 1248f3b4b4dSDmitry Karpeev PetscBool doprint,found; 1250298fd71SBarry Smith PetscViewer viewer, sviewer = NULL; 1268f3b4b4dSDmitry Karpeev PetscInt *numbering,*permutation;/* global numbering of locally-supported subdomains and the permutation from the local ordering */ 12706b43e7eSDmitry Karpeev PetscErrorCode ierr; 12806b43e7eSDmitry Karpeev 12906b43e7eSDmitry Karpeev PetscFunctionBegin; 13006b43e7eSDmitry Karpeev ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 1318f3b4b4dSDmitry Karpeev doprint = PETSC_FALSE; 1328f3b4b4dSDmitry Karpeev ierr = PetscOptionsGetBool(prefix,"-pc_gasm_print_subdomains",&doprint,NULL);CHKERRQ(ierr); 1338f3b4b4dSDmitry Karpeev if (!doprint) PetscFunctionReturn(0); 13406b43e7eSDmitry Karpeev ierr = PetscOptionsGetString(prefix,"-pc_gasm_print_subdomains",fname,PETSC_MAX_PATH_LEN,&found);CHKERRQ(ierr); 13506b43e7eSDmitry Karpeev if (!found) { ierr = PetscStrcpy(fname,"stdout");CHKERRQ(ierr); }; 136ce94432eSBarry Smith ierr = PetscViewerASCIIOpen(PetscObjectComm((PetscObject)pc),fname,&viewer);CHKERRQ(ierr); 1374bde246dSDmitry Karpeev /* 1384bde246dSDmitry Karpeev Make sure the viewer has a name. Otherwise this may cause a deadlock or other weird errors when creating a subcomm viewer: 1394bde246dSDmitry Karpeev the subcomm viewer will attempt to inherit the viewer's name, which, if not set, will be constructed collectively on the comm. 1404bde246dSDmitry Karpeev */ 1414bde246dSDmitry Karpeev ierr = PetscObjectName((PetscObject)viewer);CHKERRQ(ierr); 1424bde246dSDmitry Karpeev l = 0; 1438f3b4b4dSDmitry Karpeev ierr = PCGASMComputeGlobalSubdomainNumbering_Private(pc,&numbering,&permutation);CHKERRQ(ierr); 1448f3b4b4dSDmitry Karpeev for (count = 0; count < osm->N; ++count) { 1454bde246dSDmitry Karpeev /* Now let subdomains go one at a time in the global numbering order and print their subdomain/solver info. */ 1464bde246dSDmitry Karpeev if (l<osm->n) { 1474bde246dSDmitry Karpeev d = permutation[l]; /* d is the local number of the l-th smallest (in the global ordering) among the locally supported subdomains */ 1484bde246dSDmitry Karpeev if (numbering[d] == count) { 1493f08860eSBarry Smith ierr = PetscViewerGetSubViewer(viewer,((PetscObject)osm->ois[d])->comm, &sviewer);CHKERRQ(ierr); 1504bde246dSDmitry Karpeev ierr = PCGASMSubdomainView_Private(pc,d,sviewer);CHKERRQ(ierr); 1513f08860eSBarry Smith ierr = PetscViewerRestoreSubViewer(viewer,((PetscObject)osm->ois[d])->comm, &sviewer);CHKERRQ(ierr); 1524bde246dSDmitry Karpeev ++l; 153af538404SDmitry Karpeev } 1544bde246dSDmitry Karpeev } 155ce94432eSBarry Smith ierr = MPI_Barrier(PetscObjectComm((PetscObject)pc));CHKERRQ(ierr); 1564bde246dSDmitry Karpeev } 1578f3b4b4dSDmitry Karpeev ierr = PetscFree2(numbering,permutation);CHKERRQ(ierr); 1584bde246dSDmitry Karpeev ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 159af538404SDmitry Karpeev PetscFunctionReturn(0); 160af538404SDmitry Karpeev } 161af538404SDmitry Karpeev 162af538404SDmitry Karpeev 163af538404SDmitry Karpeev #undef __FUNCT__ 164f746d493SDmitry Karpeev #define __FUNCT__ "PCView_GASM" 165f746d493SDmitry Karpeev static PetscErrorCode PCView_GASM(PC pc,PetscViewer viewer) 166b1a0a8a3SJed Brown { 167f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 168af538404SDmitry Karpeev const char *prefix; 169b1a0a8a3SJed Brown PetscErrorCode ierr; 170af538404SDmitry Karpeev PetscMPIInt rank, size; 1718f3b4b4dSDmitry Karpeev PetscInt bsz; 17206b43e7eSDmitry Karpeev PetscBool iascii,view_subdomains=PETSC_FALSE; 173b1a0a8a3SJed Brown PetscViewer sviewer; 1748f3b4b4dSDmitry Karpeev PetscInt count, l; 1756a4f0f73SDmitry Karpeev char overlap[256] = "user-defined overlap"; 1766a4f0f73SDmitry Karpeev char gsubdomains[256] = "unknown total number of subdomains"; 17706b43e7eSDmitry Karpeev char msubdomains[256] = "unknown max number of local subdomains"; 1788f3b4b4dSDmitry Karpeev PetscInt *numbering,*permutation;/* global numbering of locally-supported subdomains and the permutation from the local ordering */ 17911aeaf0aSBarry Smith 180b1a0a8a3SJed Brown PetscFunctionBegin; 181ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size);CHKERRQ(ierr); 182ce94432eSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank);CHKERRQ(ierr); 18306b43e7eSDmitry Karpeev 18406b43e7eSDmitry Karpeev if (osm->overlap >= 0) { 18506b43e7eSDmitry Karpeev ierr = PetscSNPrintf(overlap,sizeof(overlap),"requested amount of overlap = %D",osm->overlap);CHKERRQ(ierr); 18606b43e7eSDmitry Karpeev } 1878f3b4b4dSDmitry Karpeev if (osm->N != PETSC_DETERMINE) { 1888f3b4b4dSDmitry Karpeev ierr = PetscSNPrintf(gsubdomains, sizeof(gsubdomains), "total number of subdomains = %D",osm->N);CHKERRQ(ierr); 18906b43e7eSDmitry Karpeev } 1908f3b4b4dSDmitry Karpeev if (osm->nmax != PETSC_DETERMINE) { 19106b43e7eSDmitry Karpeev ierr = PetscSNPrintf(msubdomains,sizeof(msubdomains),"max number of local subdomains = %D",osm->nmax);CHKERRQ(ierr); 19206b43e7eSDmitry Karpeev } 19306b43e7eSDmitry Karpeev 194af538404SDmitry Karpeev ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 1950298fd71SBarry Smith ierr = PetscOptionsGetBool(prefix,"-pc_gasm_view_subdomains",&view_subdomains,NULL);CHKERRQ(ierr); 19606b43e7eSDmitry Karpeev 19706b43e7eSDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 198b1a0a8a3SJed Brown if (iascii) { 19906b43e7eSDmitry Karpeev /* 20006b43e7eSDmitry Karpeev Make sure the viewer has a name. Otherwise this may cause a deadlock when creating a subcomm viewer: 20106b43e7eSDmitry Karpeev the subcomm viewer will attempt to inherit the viewer's name, which, if not set, will be constructed 20206b43e7eSDmitry Karpeev collectively on the comm. 20306b43e7eSDmitry Karpeev */ 20406b43e7eSDmitry Karpeev ierr = PetscObjectName((PetscObject)viewer);CHKERRQ(ierr); 20506b43e7eSDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer,"Generalized additive Schwarz:\n");CHKERRQ(ierr); 20606b43e7eSDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer,"Restriction/interpolation type: %s\n",PCGASMTypes[osm->type]);CHKERRQ(ierr); 20706b43e7eSDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer,"%s\n",overlap);CHKERRQ(ierr); 20806b43e7eSDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer,"%s\n",gsubdomains);CHKERRQ(ierr); 20906b43e7eSDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer,"%s\n",msubdomains);CHKERRQ(ierr); 2101575c14dSBarry Smith ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 2118f3b4b4dSDmitry Karpeev ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d|%d] number of locally-supported subdomains = %D\n",rank,size,osm->n);CHKERRQ(ierr); 212af538404SDmitry Karpeev ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 2131575c14dSBarry Smith ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 2146a4f0f73SDmitry Karpeev /* Cannot take advantage of osm->same_subdomain_solvers without a global numbering of subdomains. */ 21506b43e7eSDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer,"Subdomain solver info is as follows:\n");CHKERRQ(ierr); 216b1a0a8a3SJed Brown ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 217b1a0a8a3SJed Brown ierr = PetscViewerASCIIPrintf(viewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr); 21806b43e7eSDmitry Karpeev /* Make sure that everybody waits for the banner to be printed. */ 219ce94432eSBarry Smith ierr = MPI_Barrier(PetscObjectComm((PetscObject)viewer));CHKERRQ(ierr); 2204bde246dSDmitry Karpeev /* Now let subdomains go one at a time in the global numbering order and print their subdomain/solver info. */ 2218f3b4b4dSDmitry Karpeev ierr = PCGASMComputeGlobalSubdomainNumbering_Private(pc,&numbering,&permutation);CHKERRQ(ierr); 22206b43e7eSDmitry Karpeev l = 0; 2238f3b4b4dSDmitry Karpeev for (count = 0; count < osm->N; ++count) { 22406b43e7eSDmitry Karpeev PetscMPIInt srank, ssize; 22506b43e7eSDmitry Karpeev if (l<osm->n) { 22606b43e7eSDmitry Karpeev PetscInt d = permutation[l]; /* d is the local number of the l-th smallest (in the global ordering) among the locally supported subdomains */ 22706b43e7eSDmitry Karpeev if (numbering[d] == count) { 2286a4f0f73SDmitry Karpeev ierr = MPI_Comm_size(((PetscObject)osm->ois[d])->comm, &ssize);CHKERRQ(ierr); 2296a4f0f73SDmitry Karpeev ierr = MPI_Comm_rank(((PetscObject)osm->ois[d])->comm, &srank);CHKERRQ(ierr); 2303f08860eSBarry Smith ierr = PetscViewerGetSubViewer(viewer,((PetscObject)osm->ois[d])->comm, &sviewer);CHKERRQ(ierr); 2316a4f0f73SDmitry Karpeev ierr = ISGetLocalSize(osm->ois[d],&bsz);CHKERRQ(ierr); 2328f3b4b4dSDmitry Karpeev ierr = PetscViewerASCIISynchronizedPrintf(sviewer,"[%d|%d] (subcomm [%d|%d]) local subdomain number %D, local size = %D\n",rank,size,srank,ssize,d,bsz);CHKERRQ(ierr); 2336a4f0f73SDmitry Karpeev ierr = PetscViewerFlush(sviewer);CHKERRQ(ierr); 23406b43e7eSDmitry Karpeev if (view_subdomains) { 23506b43e7eSDmitry Karpeev ierr = PCGASMSubdomainView_Private(pc,d,sviewer);CHKERRQ(ierr); 23606b43e7eSDmitry Karpeev } 2376a4f0f73SDmitry Karpeev if (!pc->setupcalled) { 2381575c14dSBarry Smith ierr = PetscViewerASCIIPrintf(sviewer, "Solver not set up yet: PCSetUp() not yet called\n");CHKERRQ(ierr); 2392fa5cd67SKarl Rupp } else { 24006b43e7eSDmitry Karpeev ierr = KSPView(osm->ksp[d],sviewer);CHKERRQ(ierr); 2416a4f0f73SDmitry Karpeev } 24206b43e7eSDmitry Karpeev ierr = PetscViewerASCIIPrintf(sviewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr); 2437b23a99aSBarry Smith ierr = PetscViewerFlush(sviewer);CHKERRQ(ierr); 2443f08860eSBarry Smith ierr = PetscViewerRestoreSubViewer(viewer,((PetscObject)osm->ois[d])->comm, &sviewer);CHKERRQ(ierr); 24506b43e7eSDmitry Karpeev ++l; 246b1a0a8a3SJed Brown } 24706b43e7eSDmitry Karpeev } 248ce94432eSBarry Smith ierr = MPI_Barrier(PetscObjectComm((PetscObject)pc));CHKERRQ(ierr); 249b1a0a8a3SJed Brown } 2508f3b4b4dSDmitry Karpeev ierr = PetscFree2(numbering,permutation);CHKERRQ(ierr); 251b1a0a8a3SJed Brown ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 252b1a0a8a3SJed Brown } 253c5e4d11fSDmitry Karpeev ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 2541575c14dSBarry Smith ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 255b1a0a8a3SJed Brown PetscFunctionReturn(0); 256b1a0a8a3SJed Brown } 257b1a0a8a3SJed Brown 2588f3b4b4dSDmitry Karpeev PETSC_INTERN PetscErrorCode PCGASMCreateLocalSubdomains(Mat A, PetscInt nloc, IS *iis[]); 259af538404SDmitry Karpeev 260b1a0a8a3SJed Brown #undef __FUNCT__ 261f746d493SDmitry Karpeev #define __FUNCT__ "PCSetUp_GASM" 262f746d493SDmitry Karpeev static PetscErrorCode PCSetUp_GASM(PC pc) 263b1a0a8a3SJed Brown { 264f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 265b1a0a8a3SJed Brown PetscErrorCode ierr; 266b1a0a8a3SJed Brown PetscBool symset,flg; 26788389755SJed Brown PetscInt i; 268c10bc1a1SDmitry Karpeev PetscMPIInt rank, size; 269b1a0a8a3SJed Brown MatReuse scall = MAT_REUSE_MATRIX; 270b1a0a8a3SJed Brown KSP ksp; 271b1a0a8a3SJed Brown PC subpc; 272b1a0a8a3SJed Brown const char *prefix,*pprefix; 273f746d493SDmitry Karpeev Vec x,y; 2746a4f0f73SDmitry Karpeev PetscInt oni; /* Number of indices in the i-th local outer subdomain. */ 2756a4f0f73SDmitry Karpeev const PetscInt *oidxi; /* Indices from the i-th subdomain local outer subdomain. */ 2766a4f0f73SDmitry Karpeev PetscInt on; /* Number of indices in the disjoint union of local outer subdomains. */ 2776a4f0f73SDmitry Karpeev PetscInt *oidx; /* Indices in the disjoint union of local outer subdomains. */ 2786a4f0f73SDmitry Karpeev IS gois; /* Disjoint union the global indices of outer subdomains. */ 2796a4f0f73SDmitry Karpeev IS goid; /* Identity IS of the size of the disjoint union of outer subdomains. */ 280f746d493SDmitry Karpeev PetscScalar *gxarray, *gyarray; 281930d09c1SFande Kong PetscInt gostart; /* Start of locally-owned indices in the vectors -- osm->gx,osm->gy -- over the disjoint union of outer subdomains. */ 2828f3b4b4dSDmitry Karpeev PetscInt num_subdomains = 0; 2830298fd71SBarry Smith DM *subdomain_dm = NULL; 2848f3b4b4dSDmitry Karpeev char **subdomain_names = NULL; 285f771a274SFande Kong PetscInt *numbering; 2868f3b4b4dSDmitry Karpeev 287b1a0a8a3SJed Brown 288b1a0a8a3SJed Brown PetscFunctionBegin; 289ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr); 290ce94432eSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);CHKERRQ(ierr); 291b1a0a8a3SJed Brown if (!pc->setupcalled) { 292b1a0a8a3SJed Brown 293b1a0a8a3SJed Brown if (!osm->type_set) { 294b1a0a8a3SJed Brown ierr = MatIsSymmetricKnown(pc->pmat,&symset,&flg);CHKERRQ(ierr); 2952fa5cd67SKarl Rupp if (symset && flg) osm->type = PC_GASM_BASIC; 296b1a0a8a3SJed Brown } 297b1a0a8a3SJed Brown 2988f3b4b4dSDmitry Karpeev if (osm->n == PETSC_DETERMINE) { 2998f3b4b4dSDmitry Karpeev if (osm->N != PETSC_DETERMINE) { 3008f3b4b4dSDmitry Karpeev /* No local subdomains given, but the desired number of total subdomains is known, so construct them accordingly. */ 3018f3b4b4dSDmitry Karpeev ierr = PCGASMCreateSubdomains(pc->pmat,osm->N,&osm->n,&osm->iis);CHKERRQ(ierr); 3028f3b4b4dSDmitry Karpeev } else if (osm->dm_subdomains && pc->dm) { 3038f3b4b4dSDmitry Karpeev /* try pc->dm next, if allowed */ 3048f3b4b4dSDmitry Karpeev PetscInt d; 305a35b7b57SDmitry Karpeev IS *inner_subdomain_is, *outer_subdomain_is; 306a35b7b57SDmitry Karpeev ierr = DMCreateDomainDecomposition(pc->dm, &num_subdomains, &subdomain_names, &inner_subdomain_is, &outer_subdomain_is, &subdomain_dm);CHKERRQ(ierr); 307a35b7b57SDmitry Karpeev if (num_subdomains) { 308a35b7b57SDmitry Karpeev ierr = PCGASMSetSubdomains(pc, num_subdomains, inner_subdomain_is, outer_subdomain_is);CHKERRQ(ierr); 30969ca1f37SDmitry Karpeev } 310a35b7b57SDmitry Karpeev for (d = 0; d < num_subdomains; ++d) { 311a35b7b57SDmitry Karpeev if (inner_subdomain_is) {ierr = ISDestroy(&inner_subdomain_is[d]);CHKERRQ(ierr);} 312a35b7b57SDmitry Karpeev if (outer_subdomain_is) {ierr = ISDestroy(&outer_subdomain_is[d]);CHKERRQ(ierr);} 313fdc48646SDmitry Karpeev } 314a35b7b57SDmitry Karpeev ierr = PetscFree(inner_subdomain_is);CHKERRQ(ierr); 315a35b7b57SDmitry Karpeev ierr = PetscFree(outer_subdomain_is);CHKERRQ(ierr); 3168f3b4b4dSDmitry Karpeev } else { 3178f3b4b4dSDmitry Karpeev /* still no subdomains; use one per processor */ 31844fe18e5SDmitry Karpeev osm->nmax = osm->n = 1; 319ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr); 320f746d493SDmitry Karpeev osm->N = size; 3218f3b4b4dSDmitry Karpeev ierr = PCGASMCreateLocalSubdomains(pc->pmat,osm->n,&osm->iis);CHKERRQ(ierr); 322fdc48646SDmitry Karpeev } 32306b43e7eSDmitry Karpeev } 324a35b7b57SDmitry Karpeev if (!osm->iis) { 325a35b7b57SDmitry Karpeev /* 3268f3b4b4dSDmitry Karpeev osm->n was set in PCGASMSetSubdomains(), but the actual subdomains have not been supplied. 3278f3b4b4dSDmitry Karpeev We create the requisite number of local inner subdomains and then expand them into 3288f3b4b4dSDmitry Karpeev out subdomains, if necessary. 329a35b7b57SDmitry Karpeev */ 3308f3b4b4dSDmitry Karpeev ierr = PCGASMCreateLocalSubdomains(pc->pmat,osm->n,&osm->iis);CHKERRQ(ierr); 331a35b7b57SDmitry Karpeev } 3328f3b4b4dSDmitry Karpeev if (!osm->ois) { 3338f3b4b4dSDmitry Karpeev /* 3348f3b4b4dSDmitry Karpeev Initially make outer subdomains the same as inner subdomains. If nonzero additional overlap 3358f3b4b4dSDmitry Karpeev has been requested, copy the inner subdomains over so they can be modified. 3368f3b4b4dSDmitry Karpeev */ 3378f3b4b4dSDmitry Karpeev ierr = PetscMalloc1(osm->n,&osm->ois);CHKERRQ(ierr); 3388f3b4b4dSDmitry Karpeev for (i=0; i<osm->n; ++i) { 3398f3b4b4dSDmitry Karpeev if (osm->overlap > 0) { /* With positive overlap, osm->iis[i] will be modified */ 3408f3b4b4dSDmitry Karpeev ierr = ISDuplicate(osm->iis[i],(osm->ois)+i);CHKERRQ(ierr); 3418f3b4b4dSDmitry Karpeev ierr = ISCopy(osm->iis[i],osm->ois[i]);CHKERRQ(ierr); 3428f3b4b4dSDmitry Karpeev } else { 3438f3b4b4dSDmitry Karpeev ierr = PetscObjectReference((PetscObject)((osm->iis)[i]));CHKERRQ(ierr); 3448f3b4b4dSDmitry Karpeev osm->ois[i] = osm->iis[i]; 3458f3b4b4dSDmitry Karpeev } 3468f3b4b4dSDmitry Karpeev } 3478f3b4b4dSDmitry Karpeev if (osm->overlap > 0) { 3488f3b4b4dSDmitry Karpeev /* Extend the "overlapping" regions by a number of steps */ 349d21f2a47SFande Kong ierr = MatIncreaseOverlapSplit(pc->pmat,osm->n,osm->ois,osm->overlap);CHKERRQ(ierr); 3508f3b4b4dSDmitry Karpeev } 351b1a0a8a3SJed Brown } 3526a4f0f73SDmitry Karpeev 3538f3b4b4dSDmitry Karpeev /* Now the subdomains are defined. Determine their global and max local numbers, if necessary. */ 3548f3b4b4dSDmitry Karpeev if (osm->nmax == PETSC_DETERMINE) { 3558f3b4b4dSDmitry Karpeev PetscMPIInt inwork,outwork; 3568f3b4b4dSDmitry Karpeev /* determine global number of subdomains and the max number of local subdomains */ 3578f3b4b4dSDmitry Karpeev inwork = osm->n; 358*b2566f29SBarry Smith ierr = MPIU_Allreduce(&inwork,&outwork,1,MPI_INT,MPI_MAX,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr); 3598f3b4b4dSDmitry Karpeev osm->nmax = outwork; 3608f3b4b4dSDmitry Karpeev } 3618f3b4b4dSDmitry Karpeev if (osm->N == PETSC_DETERMINE) { 3628f3b4b4dSDmitry Karpeev /* Determine the number of globally-distinct subdomains and compute a global numbering for them. */ 3638f3b4b4dSDmitry Karpeev ierr = PetscObjectsListGetGlobalNumbering(PetscObjectComm((PetscObject)pc),osm->n,(PetscObject*)osm->ois,&osm->N,NULL);CHKERRQ(ierr); 3648f3b4b4dSDmitry Karpeev } 3658f3b4b4dSDmitry Karpeev 366b1a0a8a3SJed Brown 367b1a0a8a3SJed Brown if (osm->sort_indices) { 368f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 3696a4f0f73SDmitry Karpeev ierr = ISSort(osm->ois[i]);CHKERRQ(ierr); 3706a4f0f73SDmitry Karpeev ierr = ISSort(osm->iis[i]);CHKERRQ(ierr); 371b1a0a8a3SJed Brown } 372b1a0a8a3SJed Brown } 3738f3b4b4dSDmitry Karpeev ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 3748f3b4b4dSDmitry Karpeev ierr = PCGASMPrintSubdomains(pc);CHKERRQ(ierr); 3758f3b4b4dSDmitry Karpeev 3766a4f0f73SDmitry Karpeev /* 3776a4f0f73SDmitry Karpeev Merge the ISs, create merged vectors and restrictions. 3786a4f0f73SDmitry Karpeev */ 3796a4f0f73SDmitry Karpeev /* Merge outer subdomain ISs and construct a restriction onto the disjoint union of local outer subdomains. */ 3806a4f0f73SDmitry Karpeev on = 0; 381f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 3826a4f0f73SDmitry Karpeev ierr = ISGetLocalSize(osm->ois[i],&oni);CHKERRQ(ierr); 3836a4f0f73SDmitry Karpeev on += oni; 384f746d493SDmitry Karpeev } 385785e854fSJed Brown ierr = PetscMalloc1(on, &oidx);CHKERRQ(ierr); 3866a4f0f73SDmitry Karpeev on = 0; 387930d09c1SFande Kong /* Merge local indices together */ 388f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 3896a4f0f73SDmitry Karpeev ierr = ISGetLocalSize(osm->ois[i],&oni);CHKERRQ(ierr); 3906a4f0f73SDmitry Karpeev ierr = ISGetIndices(osm->ois[i],&oidxi);CHKERRQ(ierr); 3916a4f0f73SDmitry Karpeev ierr = PetscMemcpy(oidx+on,oidxi,sizeof(PetscInt)*oni);CHKERRQ(ierr); 3926a4f0f73SDmitry Karpeev ierr = ISRestoreIndices(osm->ois[i],&oidxi);CHKERRQ(ierr); 3936a4f0f73SDmitry Karpeev on += oni; 394f746d493SDmitry Karpeev } 3956a4f0f73SDmitry Karpeev ierr = ISCreateGeneral(((PetscObject)(pc))->comm,on,oidx,PETSC_OWN_POINTER,&gois);CHKERRQ(ierr); 3962a7a6963SBarry Smith ierr = MatCreateVecs(pc->pmat,&x,&y);CHKERRQ(ierr); 397ce94432eSBarry Smith ierr = VecCreateMPI(PetscObjectComm((PetscObject)pc),on,PETSC_DECIDE,&osm->gx);CHKERRQ(ierr); 398f746d493SDmitry Karpeev ierr = VecDuplicate(osm->gx,&osm->gy);CHKERRQ(ierr); 399930d09c1SFande Kong ierr = VecGetOwnershipRange(osm->gx, &gostart, NULL);CHKERRQ(ierr); 400930d09c1SFande Kong ierr = ISCreateStride(PetscObjectComm((PetscObject)pc),on,gostart,1, &goid);CHKERRQ(ierr); 401930d09c1SFande Kong /* gois might indices not on local */ 4026a4f0f73SDmitry Karpeev ierr = VecScatterCreate(x,gois,osm->gx,goid, &(osm->gorestriction));CHKERRQ(ierr); 403f771a274SFande Kong ierr = PetscMalloc1(osm->n,&numbering);CHKERRQ(ierr); 404f771a274SFande Kong ierr = PetscObjectsListGetGlobalNumbering(PetscObjectComm((PetscObject)pc),osm->n,(PetscObject*)osm->ois,NULL,numbering);CHKERRQ(ierr); 4056a4f0f73SDmitry Karpeev ierr = VecDestroy(&x);CHKERRQ(ierr); 4067105d80bSDmitry Karpeev ierr = ISDestroy(&gois);CHKERRQ(ierr); 4072fa5cd67SKarl Rupp 4086a4f0f73SDmitry Karpeev /* Merge inner subdomain ISs and construct a restriction onto the disjoint union of local inner subdomains. */ 4092fa5cd67SKarl Rupp { 4102fa5cd67SKarl Rupp PetscInt ini; /* Number of indices the i-th a local inner subdomain. */ 4116a4f0f73SDmitry Karpeev PetscInt in; /* Number of indices in the disjoint uniont of local inner subdomains. */ 4126a4f0f73SDmitry Karpeev PetscInt *iidx; /* Global indices in the merged local inner subdomain. */ 4136a4f0f73SDmitry Karpeev PetscInt *ioidx; /* Global indices of the disjoint union of inner subdomains within the disjoint union of outer subdomains. */ 4146a4f0f73SDmitry Karpeev IS giis; /* IS for the disjoint union of inner subdomains. */ 4156a4f0f73SDmitry Karpeev IS giois; /* IS for the disjoint union of inner subdomains within the disjoint union of outer subdomains. */ 416f771a274SFande Kong PetscScalar *array; 417f771a274SFande Kong const PetscInt *indices; 418f771a274SFande Kong PetscInt k; 4196a4f0f73SDmitry Karpeev on = 0; 420f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 4216a4f0f73SDmitry Karpeev ierr = ISGetLocalSize(osm->ois[i],&oni);CHKERRQ(ierr); 4226a4f0f73SDmitry Karpeev on += oni; 423f746d493SDmitry Karpeev } 424f771a274SFande Kong ierr = PetscMalloc1(on, &iidx);CHKERRQ(ierr); 425f771a274SFande Kong ierr = PetscMalloc1(on, &ioidx);CHKERRQ(ierr); 426f771a274SFande Kong ierr = VecGetArray(y,&array);CHKERRQ(ierr); 427e12b4729SFande Kong /* set communicator id to determine where overlap is */ 428f771a274SFande Kong in = 0; 429f771a274SFande Kong for (i=0; i<osm->n; i++) { 430f771a274SFande Kong ierr = ISGetLocalSize(osm->iis[i],&ini);CHKERRQ(ierr); 431f771a274SFande Kong for (k = 0; k < ini; ++k){ 432f771a274SFande Kong array[in+k] = numbering[i]; 433f771a274SFande Kong } 434f771a274SFande Kong in += ini; 435f771a274SFande Kong } 436f771a274SFande Kong ierr = VecRestoreArray(y,&array);CHKERRQ(ierr); 437f771a274SFande Kong ierr = VecScatterBegin(osm->gorestriction,y,osm->gy,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 438f771a274SFande Kong ierr = VecScatterEnd(osm->gorestriction,y,osm->gy,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 439f771a274SFande Kong ierr = VecGetOwnershipRange(osm->gy,&gostart, NULL);CHKERRQ(ierr); 440f771a274SFande Kong ierr = VecGetArray(osm->gy,&array);CHKERRQ(ierr); 441f771a274SFande Kong on = 0; 442f771a274SFande Kong in = 0; 443f771a274SFande Kong for (i=0; i<osm->n; i++) { 444f771a274SFande Kong ierr = ISGetLocalSize(osm->ois[i],&oni);CHKERRQ(ierr); 445f771a274SFande Kong ierr = ISGetIndices(osm->ois[i],&indices);CHKERRQ(ierr); 446f771a274SFande Kong for (k=0; k<oni; k++) { 447e12b4729SFande Kong /* skip overlapping indices to get inner domain */ 44843081b6cSDmitry Karpeev if(PetscRealPart(array[on+k]) != numbering[i]) continue; 449f771a274SFande Kong iidx[in] = indices[k]; 450f771a274SFande Kong ioidx[in++] = gostart+on+k; 451f771a274SFande Kong } 452f771a274SFande Kong ierr = ISRestoreIndices(osm->ois[i], &indices);CHKERRQ(ierr); 453f771a274SFande Kong on += oni; 454f771a274SFande Kong } 455f771a274SFande Kong ierr = VecRestoreArray(osm->gy,&array);CHKERRQ(ierr); 456ce94432eSBarry Smith ierr = ISCreateGeneral(PetscObjectComm((PetscObject)pc),in,iidx,PETSC_OWN_POINTER,&giis);CHKERRQ(ierr); 457ce94432eSBarry Smith ierr = ISCreateGeneral(PetscObjectComm((PetscObject)pc),in,ioidx,PETSC_OWN_POINTER,&giois);CHKERRQ(ierr); 4586a4f0f73SDmitry Karpeev ierr = VecScatterCreate(y,giis,osm->gy,giois,&osm->girestriction);CHKERRQ(ierr); 4596a4f0f73SDmitry Karpeev ierr = VecDestroy(&y);CHKERRQ(ierr); 4606a4f0f73SDmitry Karpeev ierr = ISDestroy(&giis);CHKERRQ(ierr); 4616a4f0f73SDmitry Karpeev ierr = ISDestroy(&giois);CHKERRQ(ierr); 462b1a0a8a3SJed Brown } 4636a4f0f73SDmitry Karpeev ierr = ISDestroy(&goid);CHKERRQ(ierr); 464f771a274SFande Kong ierr = PetscFree(numbering);CHKERRQ(ierr); 4652fa5cd67SKarl Rupp 466f746d493SDmitry Karpeev /* Create the subdomain work vectors. */ 467785e854fSJed Brown ierr = PetscMalloc1(osm->n,&osm->x);CHKERRQ(ierr); 468785e854fSJed Brown ierr = PetscMalloc1(osm->n,&osm->y);CHKERRQ(ierr); 469f746d493SDmitry Karpeev ierr = VecGetArray(osm->gx, &gxarray);CHKERRQ(ierr); 470f746d493SDmitry Karpeev ierr = VecGetArray(osm->gy, &gyarray);CHKERRQ(ierr); 4716a4f0f73SDmitry Karpeev for (i=0, on=0; i<osm->n; ++i, on += oni) { 4726a4f0f73SDmitry Karpeev PetscInt oNi; 4736a4f0f73SDmitry Karpeev ierr = ISGetLocalSize(osm->ois[i],&oni);CHKERRQ(ierr); 474930d09c1SFande Kong /* on a sub communicator */ 4756a4f0f73SDmitry Karpeev ierr = ISGetSize(osm->ois[i],&oNi);CHKERRQ(ierr); 4766a4f0f73SDmitry Karpeev ierr = VecCreateMPIWithArray(((PetscObject)(osm->ois[i]))->comm,1,oni,oNi,gxarray+on,&osm->x[i]);CHKERRQ(ierr); 4776a4f0f73SDmitry Karpeev ierr = VecCreateMPIWithArray(((PetscObject)(osm->ois[i]))->comm,1,oni,oNi,gyarray+on,&osm->y[i]);CHKERRQ(ierr); 478b1a0a8a3SJed Brown } 479f746d493SDmitry Karpeev ierr = VecRestoreArray(osm->gx, &gxarray);CHKERRQ(ierr); 480f746d493SDmitry Karpeev ierr = VecRestoreArray(osm->gy, &gyarray);CHKERRQ(ierr); 4818f3b4b4dSDmitry Karpeev /* Create the subdomain solvers */ 482785e854fSJed Brown ierr = PetscMalloc1(osm->n,&osm->ksp);CHKERRQ(ierr); 4838f3b4b4dSDmitry Karpeev for (i=0; i<osm->n; i++) { 4848f3b4b4dSDmitry Karpeev char subprefix[PETSC_MAX_PATH_LEN+1]; 4856a4f0f73SDmitry Karpeev ierr = KSPCreate(((PetscObject)(osm->ois[i]))->comm,&ksp);CHKERRQ(ierr); 4863bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ksp);CHKERRQ(ierr); 487b1a0a8a3SJed Brown ierr = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,1);CHKERRQ(ierr); 488b1a0a8a3SJed Brown ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr); 4898f3b4b4dSDmitry Karpeev ierr = KSPGetPC(ksp,&subpc);CHKERRQ(ierr); /* Why do we need this here? */ 4908f3b4b4dSDmitry Karpeev if (subdomain_dm) { 4918f3b4b4dSDmitry Karpeev ierr = KSPSetDM(ksp,subdomain_dm[i]);CHKERRQ(ierr); 4928f3b4b4dSDmitry Karpeev ierr = DMDestroy(subdomain_dm+i);CHKERRQ(ierr); 4938f3b4b4dSDmitry Karpeev } 494b1a0a8a3SJed Brown ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 495b1a0a8a3SJed Brown ierr = KSPSetOptionsPrefix(ksp,prefix);CHKERRQ(ierr); 4968f3b4b4dSDmitry Karpeev if (subdomain_names && subdomain_names[i]) { 4978f3b4b4dSDmitry Karpeev ierr = PetscSNPrintf(subprefix,PETSC_MAX_PATH_LEN,"sub_%s_",subdomain_names[i]);CHKERRQ(ierr); 4988f3b4b4dSDmitry Karpeev ierr = KSPAppendOptionsPrefix(ksp,subprefix);CHKERRQ(ierr); 4998f3b4b4dSDmitry Karpeev ierr = PetscFree(subdomain_names[i]);CHKERRQ(ierr); 5008f3b4b4dSDmitry Karpeev } 501b1a0a8a3SJed Brown ierr = KSPAppendOptionsPrefix(ksp,"sub_");CHKERRQ(ierr); 502b1a0a8a3SJed Brown osm->ksp[i] = ksp; 503b1a0a8a3SJed Brown } 5048f3b4b4dSDmitry Karpeev ierr = PetscFree(subdomain_dm);CHKERRQ(ierr); 5058f3b4b4dSDmitry Karpeev ierr = PetscFree(subdomain_names);CHKERRQ(ierr); 506b1a0a8a3SJed Brown scall = MAT_INITIAL_MATRIX; 507b1a0a8a3SJed Brown 5088f3b4b4dSDmitry Karpeev } else { /* if (pc->setupcalled) */ 509b1a0a8a3SJed Brown /* 5108f3b4b4dSDmitry Karpeev Destroy the submatrices from the previous iteration 511b1a0a8a3SJed Brown */ 512b1a0a8a3SJed Brown if (pc->flag == DIFFERENT_NONZERO_PATTERN) { 513f746d493SDmitry Karpeev ierr = MatDestroyMatrices(osm->n,&osm->pmat);CHKERRQ(ierr); 514b1a0a8a3SJed Brown scall = MAT_INITIAL_MATRIX; 515b1a0a8a3SJed Brown } 516b1a0a8a3SJed Brown } 517b1a0a8a3SJed Brown 518b1a0a8a3SJed Brown /* 519f746d493SDmitry Karpeev Extract out the submatrices. 520b1a0a8a3SJed Brown */ 52181a5c4aaSDmitry Karpeev if (size > 1) { 5228f3b4b4dSDmitry Karpeev ierr = MatGetSubMatricesMPI(pc->pmat,osm->n,osm->ois,osm->ois,scall,&osm->pmat);CHKERRQ(ierr); 5232fa5cd67SKarl Rupp } else { 5246a4f0f73SDmitry Karpeev ierr = MatGetSubMatrices(pc->pmat,osm->n,osm->ois,osm->ois,scall,&osm->pmat);CHKERRQ(ierr); 52581a5c4aaSDmitry Karpeev } 526b1a0a8a3SJed Brown if (scall == MAT_INITIAL_MATRIX) { 527b1a0a8a3SJed Brown ierr = PetscObjectGetOptionsPrefix((PetscObject)pc->pmat,&pprefix);CHKERRQ(ierr); 528f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 5293bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)osm->pmat[i]);CHKERRQ(ierr); 530b1a0a8a3SJed Brown ierr = PetscObjectSetOptionsPrefix((PetscObject)osm->pmat[i],pprefix);CHKERRQ(ierr); 531b1a0a8a3SJed Brown } 532b1a0a8a3SJed Brown } 533b1a0a8a3SJed Brown 534b1a0a8a3SJed Brown /* Return control to the user so that the submatrices can be modified (e.g., to apply 535b1a0a8a3SJed Brown different boundary conditions for the submatrices than for the global problem) */ 5366a4f0f73SDmitry Karpeev ierr = PCModifySubMatrices(pc,osm->n,osm->ois,osm->ois,osm->pmat,pc->modifysubmatricesP);CHKERRQ(ierr); 537b1a0a8a3SJed Brown 538b1a0a8a3SJed Brown /* 5396a4f0f73SDmitry Karpeev Loop over submatrices putting them into local ksps 540b1a0a8a3SJed Brown */ 541f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 54223ee1639SBarry Smith ierr = KSPSetOperators(osm->ksp[i],osm->pmat[i],osm->pmat[i]);CHKERRQ(ierr); 543b1a0a8a3SJed Brown if (!pc->setupcalled) { 544b1a0a8a3SJed Brown ierr = KSPSetFromOptions(osm->ksp[i]);CHKERRQ(ierr); 545b1a0a8a3SJed Brown } 546b1a0a8a3SJed Brown } 547b1a0a8a3SJed Brown PetscFunctionReturn(0); 548b1a0a8a3SJed Brown } 549b1a0a8a3SJed Brown 550b1a0a8a3SJed Brown #undef __FUNCT__ 551f746d493SDmitry Karpeev #define __FUNCT__ "PCSetUpOnBlocks_GASM" 552f746d493SDmitry Karpeev static PetscErrorCode PCSetUpOnBlocks_GASM(PC pc) 553b1a0a8a3SJed Brown { 554f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 555b1a0a8a3SJed Brown PetscErrorCode ierr; 556b1a0a8a3SJed Brown PetscInt i; 557b1a0a8a3SJed Brown 558b1a0a8a3SJed Brown PetscFunctionBegin; 559f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 560b1a0a8a3SJed Brown ierr = KSPSetUp(osm->ksp[i]);CHKERRQ(ierr); 561b1a0a8a3SJed Brown } 562b1a0a8a3SJed Brown PetscFunctionReturn(0); 563b1a0a8a3SJed Brown } 564b1a0a8a3SJed Brown 565b1a0a8a3SJed Brown #undef __FUNCT__ 566f746d493SDmitry Karpeev #define __FUNCT__ "PCApply_GASM" 567f746d493SDmitry Karpeev static PetscErrorCode PCApply_GASM(PC pc,Vec x,Vec y) 568b1a0a8a3SJed Brown { 569f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 570b1a0a8a3SJed Brown PetscErrorCode ierr; 571f746d493SDmitry Karpeev PetscInt i; 572b1a0a8a3SJed Brown ScatterMode forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE; 573b1a0a8a3SJed Brown 574b1a0a8a3SJed Brown PetscFunctionBegin; 575b1a0a8a3SJed Brown /* 5766a4f0f73SDmitry Karpeev Support for limiting the restriction or interpolation only to the inner 577b1a0a8a3SJed Brown subdomain values (leaving the other values 0). 578b1a0a8a3SJed Brown */ 579f746d493SDmitry Karpeev if (!(osm->type & PC_GASM_RESTRICT)) { 580b1a0a8a3SJed Brown /* have to zero the work RHS since scatter may leave some slots empty */ 581f746d493SDmitry Karpeev ierr = VecZeroEntries(osm->gx);CHKERRQ(ierr); 5826a4f0f73SDmitry Karpeev ierr = VecScatterBegin(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 5832fa5cd67SKarl Rupp } else { 5846a4f0f73SDmitry Karpeev ierr = VecScatterBegin(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 585b1a0a8a3SJed Brown } 5866a4f0f73SDmitry Karpeev ierr = VecZeroEntries(osm->gy);CHKERRQ(ierr); 5876a4f0f73SDmitry Karpeev if (!(osm->type & PC_GASM_RESTRICT)) { 5886a4f0f73SDmitry Karpeev ierr = VecScatterEnd(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 5892fa5cd67SKarl Rupp } else { 5906a4f0f73SDmitry Karpeev ierr = VecScatterEnd(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 5916a4f0f73SDmitry Karpeev } 59286b96d36SDmitry Karpeev /* do the subdomain solves */ 59386b96d36SDmitry Karpeev for (i=0; i<osm->n; ++i) { 594b1a0a8a3SJed Brown ierr = KSPSolve(osm->ksp[i],osm->x[i],osm->y[i]);CHKERRQ(ierr); 595b1a0a8a3SJed Brown } 596930d09c1SFande Kong /* Do we need to zero y ?? */ 5976a4f0f73SDmitry Karpeev ierr = VecZeroEntries(y);CHKERRQ(ierr); 5986a4f0f73SDmitry Karpeev if (!(osm->type & PC_GASM_INTERPOLATE)) { 5996a4f0f73SDmitry Karpeev ierr = VecScatterBegin(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); 6006a4f0f73SDmitry Karpeev ierr = VecScatterEnd(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); PetscFunctionReturn(0); 6012fa5cd67SKarl Rupp } else { 6026a4f0f73SDmitry Karpeev ierr = VecScatterBegin(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); 6036a4f0f73SDmitry Karpeev ierr = VecScatterEnd(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); PetscFunctionReturn(0); 6046a4f0f73SDmitry Karpeev } 605b1a0a8a3SJed Brown } 606b1a0a8a3SJed Brown 607b1a0a8a3SJed Brown #undef __FUNCT__ 608f746d493SDmitry Karpeev #define __FUNCT__ "PCApplyTranspose_GASM" 609f746d493SDmitry Karpeev static PetscErrorCode PCApplyTranspose_GASM(PC pc,Vec x,Vec y) 610b1a0a8a3SJed Brown { 611f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 612b1a0a8a3SJed Brown PetscErrorCode ierr; 613f746d493SDmitry Karpeev PetscInt i; 614b1a0a8a3SJed Brown ScatterMode forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE; 615b1a0a8a3SJed Brown 616b1a0a8a3SJed Brown PetscFunctionBegin; 617b1a0a8a3SJed Brown /* 618b1a0a8a3SJed Brown Support for limiting the restriction or interpolation to only local 619b1a0a8a3SJed Brown subdomain values (leaving the other values 0). 620b1a0a8a3SJed Brown 621f746d493SDmitry Karpeev Note: these are reversed from the PCApply_GASM() because we are applying the 622b1a0a8a3SJed Brown transpose of the three terms 623b1a0a8a3SJed Brown */ 624f746d493SDmitry Karpeev if (!(osm->type & PC_GASM_INTERPOLATE)) { 625b1a0a8a3SJed Brown /* have to zero the work RHS since scatter may leave some slots empty */ 626f746d493SDmitry Karpeev ierr = VecZeroEntries(osm->gx);CHKERRQ(ierr); 6276a4f0f73SDmitry Karpeev ierr = VecScatterBegin(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 6282fa5cd67SKarl Rupp } else { 6296a4f0f73SDmitry Karpeev ierr = VecScatterBegin(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 630b1a0a8a3SJed Brown } 6316a4f0f73SDmitry Karpeev ierr = VecZeroEntries(osm->gy);CHKERRQ(ierr); 6326a4f0f73SDmitry Karpeev if (!(osm->type & PC_GASM_INTERPOLATE)) { 6336a4f0f73SDmitry Karpeev ierr = VecScatterEnd(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 6342fa5cd67SKarl Rupp } else { 6356a4f0f73SDmitry Karpeev ierr = VecScatterEnd(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 6366a4f0f73SDmitry Karpeev } 637b1a0a8a3SJed Brown /* do the local solves */ 638ab3e923fSDmitry Karpeev for (i=0; i<osm->n; ++i) { /* Note that the solves are local, so we can go to osm->n, rather than osm->nmax. */ 639b1a0a8a3SJed Brown ierr = KSPSolveTranspose(osm->ksp[i],osm->x[i],osm->y[i]);CHKERRQ(ierr); 640b1a0a8a3SJed Brown } 6416a4f0f73SDmitry Karpeev ierr = VecZeroEntries(y);CHKERRQ(ierr); 6426a4f0f73SDmitry Karpeev if (!(osm->type & PC_GASM_RESTRICT)) { 6436a4f0f73SDmitry Karpeev ierr = VecScatterBegin(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); 6446a4f0f73SDmitry Karpeev ierr = VecScatterEnd(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); 6452fa5cd67SKarl Rupp } else { 6466a4f0f73SDmitry Karpeev ierr = VecScatterBegin(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); 6476a4f0f73SDmitry Karpeev ierr = VecScatterEnd(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); 6486a4f0f73SDmitry Karpeev } 649b1a0a8a3SJed Brown PetscFunctionReturn(0); 650b1a0a8a3SJed Brown } 651b1a0a8a3SJed Brown 652b1a0a8a3SJed Brown #undef __FUNCT__ 653a06653b4SBarry Smith #define __FUNCT__ "PCReset_GASM" 654a06653b4SBarry Smith static PetscErrorCode PCReset_GASM(PC pc) 655b1a0a8a3SJed Brown { 656f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 657b1a0a8a3SJed Brown PetscErrorCode ierr; 658b1a0a8a3SJed Brown PetscInt i; 659b1a0a8a3SJed Brown 660b1a0a8a3SJed Brown PetscFunctionBegin; 661b1a0a8a3SJed Brown if (osm->ksp) { 662f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 663a06653b4SBarry Smith ierr = KSPReset(osm->ksp[i]);CHKERRQ(ierr); 664b1a0a8a3SJed Brown } 665b1a0a8a3SJed Brown } 666b1a0a8a3SJed Brown if (osm->pmat) { 667f746d493SDmitry Karpeev if (osm->n > 0) { 668ab3e923fSDmitry Karpeev ierr = MatDestroyMatrices(osm->n,&osm->pmat);CHKERRQ(ierr); 669b1a0a8a3SJed Brown } 670b1a0a8a3SJed Brown } 671d34fcf5fSBarry Smith if (osm->x) { 672f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 673fcfd50ebSBarry Smith ierr = VecDestroy(&osm->x[i]);CHKERRQ(ierr); 674fcfd50ebSBarry Smith ierr = VecDestroy(&osm->y[i]);CHKERRQ(ierr); 675b1a0a8a3SJed Brown } 676d34fcf5fSBarry Smith } 6776bf464f9SBarry Smith ierr = VecDestroy(&osm->gx);CHKERRQ(ierr); 6786bf464f9SBarry Smith ierr = VecDestroy(&osm->gy);CHKERRQ(ierr); 679ab3e923fSDmitry Karpeev 6806a4f0f73SDmitry Karpeev ierr = VecScatterDestroy(&osm->gorestriction);CHKERRQ(ierr); 6816a4f0f73SDmitry Karpeev ierr = VecScatterDestroy(&osm->girestriction);CHKERRQ(ierr); 6828f3b4b4dSDmitry Karpeev if (!osm->user_subdomains) { 6832c112581SDmitry Karpeev ierr = PCGASMDestroySubdomains(osm->n,&osm->ois,&osm->iis);CHKERRQ(ierr); 6848f3b4b4dSDmitry Karpeev osm->N = PETSC_DETERMINE; 6858f3b4b4dSDmitry Karpeev osm->nmax = PETSC_DETERMINE; 6868f3b4b4dSDmitry Karpeev } 687a06653b4SBarry Smith PetscFunctionReturn(0); 688a06653b4SBarry Smith } 689a06653b4SBarry Smith 690a06653b4SBarry Smith #undef __FUNCT__ 691a06653b4SBarry Smith #define __FUNCT__ "PCDestroy_GASM" 692a06653b4SBarry Smith static PetscErrorCode PCDestroy_GASM(PC pc) 693a06653b4SBarry Smith { 694a06653b4SBarry Smith PC_GASM *osm = (PC_GASM*)pc->data; 695a06653b4SBarry Smith PetscErrorCode ierr; 696a06653b4SBarry Smith PetscInt i; 697a06653b4SBarry Smith 698a06653b4SBarry Smith PetscFunctionBegin; 699a06653b4SBarry Smith ierr = PCReset_GASM(pc);CHKERRQ(ierr); 700135757f6SDmitry Karpeev 701135757f6SDmitry Karpeev /* PCReset will not destroy subdomains, if user_subdomains is true. */ 702135757f6SDmitry Karpeev ierr = PCGASMDestroySubdomains(osm->n,&osm->ois,&osm->iis);CHKERRQ(ierr); 703135757f6SDmitry Karpeev 704a06653b4SBarry Smith if (osm->ksp) { 705a06653b4SBarry Smith for (i=0; i<osm->n; i++) { 7066bf464f9SBarry Smith ierr = KSPDestroy(&osm->ksp[i]);CHKERRQ(ierr); 707a06653b4SBarry Smith } 708a06653b4SBarry Smith ierr = PetscFree(osm->ksp);CHKERRQ(ierr); 709a06653b4SBarry Smith } 710a06653b4SBarry Smith ierr = PetscFree(osm->x);CHKERRQ(ierr); 711a06653b4SBarry Smith ierr = PetscFree(osm->y);CHKERRQ(ierr); 712c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 713b1a0a8a3SJed Brown PetscFunctionReturn(0); 714b1a0a8a3SJed Brown } 715b1a0a8a3SJed Brown 716b1a0a8a3SJed Brown #undef __FUNCT__ 717f746d493SDmitry Karpeev #define __FUNCT__ "PCSetFromOptions_GASM" 7188c34d3f5SBarry Smith static PetscErrorCode PCSetFromOptions_GASM(PetscOptions *PetscOptionsObject,PC pc) 719a6dfd86eSKarl Rupp { 720f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 721b1a0a8a3SJed Brown PetscErrorCode ierr; 722b1a0a8a3SJed Brown PetscInt blocks,ovl; 723b1a0a8a3SJed Brown PetscBool symset,flg; 724f746d493SDmitry Karpeev PCGASMType gasmtype; 725b1a0a8a3SJed Brown 726b1a0a8a3SJed Brown PetscFunctionBegin; 727b1a0a8a3SJed Brown /* set the type to symmetric if matrix is symmetric */ 728b1a0a8a3SJed Brown if (!osm->type_set && pc->pmat) { 729b1a0a8a3SJed Brown ierr = MatIsSymmetricKnown(pc->pmat,&symset,&flg);CHKERRQ(ierr); 7302fa5cd67SKarl Rupp if (symset && flg) osm->type = PC_GASM_BASIC; 731b1a0a8a3SJed Brown } 732e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"Generalized additive Schwarz options");CHKERRQ(ierr); 7338f3b4b4dSDmitry Karpeev ierr = PetscOptionsBool("-pc_gasm_use_dm_subdomains","If subdomains aren't set, use DMCreateDomainDecomposition() to define subdomains.","PCGASMSetUseDMSubdomains",osm->dm_subdomains,&osm->dm_subdomains,&flg);CHKERRQ(ierr); 7348f3b4b4dSDmitry Karpeev ierr = PetscOptionsInt("-pc_gasm_total_subdomains","Total number of subdomains across communicator","PCGASMSetTotalSubdomains",osm->N,&blocks,&flg);CHKERRQ(ierr); 73565db9045SDmitry Karpeev if (flg) { 7368f3b4b4dSDmitry Karpeev ierr = PCGASMSetTotalSubdomains(pc,blocks);CHKERRQ(ierr); 73765db9045SDmitry Karpeev } 73806b43e7eSDmitry Karpeev ierr = PetscOptionsInt("-pc_gasm_overlap","Number of overlapping degrees of freedom","PCGASMSetOverlap",osm->overlap,&ovl,&flg);CHKERRQ(ierr); 73965db9045SDmitry Karpeev if (flg) { 74065db9045SDmitry Karpeev ierr = PCGASMSetOverlap(pc,ovl);CHKERRQ(ierr); 741d709ea83SDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 74265db9045SDmitry Karpeev } 743b1a0a8a3SJed Brown flg = PETSC_FALSE; 744f746d493SDmitry Karpeev ierr = PetscOptionsEnum("-pc_gasm_type","Type of restriction/extension","PCGASMSetType",PCGASMTypes,(PetscEnum)osm->type,(PetscEnum*)&gasmtype,&flg);CHKERRQ(ierr); 745f746d493SDmitry Karpeev if (flg) {ierr = PCGASMSetType(pc,gasmtype);CHKERRQ(ierr);} 746b1a0a8a3SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 747b1a0a8a3SJed Brown PetscFunctionReturn(0); 748b1a0a8a3SJed Brown } 749b1a0a8a3SJed Brown 750b1a0a8a3SJed Brown /*------------------------------------------------------------------------------------*/ 751b1a0a8a3SJed Brown 752b1a0a8a3SJed Brown #undef __FUNCT__ 7538f3b4b4dSDmitry Karpeev #define __FUNCT__ "PCGASMSetTotalSubdomains" 7548f3b4b4dSDmitry Karpeev /*@ 7558f3b4b4dSDmitry Karpeev PCGASMSetTotalSubdomains - sets the total number of subdomains to use across the 7568f3b4b4dSDmitry Karpeev communicator. 7578f3b4b4dSDmitry Karpeev Logically collective on pc 7588f3b4b4dSDmitry Karpeev 7598f3b4b4dSDmitry Karpeev Input Parameters: 7608f3b4b4dSDmitry Karpeev + pc - the preconditioner 7618f3b4b4dSDmitry Karpeev - N - total number of subdomains 7628f3b4b4dSDmitry Karpeev 7638f3b4b4dSDmitry Karpeev 7648f3b4b4dSDmitry Karpeev Level: beginner 7658f3b4b4dSDmitry Karpeev 7668f3b4b4dSDmitry Karpeev .keywords: PC, ASM, DM, set, subdomains, additive Schwarz 7678f3b4b4dSDmitry Karpeev 7688f3b4b4dSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMSetOverlap() 7698f3b4b4dSDmitry Karpeev PCGASMCreateSubdomains2D() 7708f3b4b4dSDmitry Karpeev @*/ 7718f3b4b4dSDmitry Karpeev PetscErrorCode PCGASMSetTotalSubdomains(PC pc,PetscInt N) 7728f3b4b4dSDmitry Karpeev { 7738f3b4b4dSDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 7748f3b4b4dSDmitry Karpeev PetscMPIInt size,rank; 7758f3b4b4dSDmitry Karpeev PetscErrorCode ierr; 7768f3b4b4dSDmitry Karpeev 7778f3b4b4dSDmitry Karpeev PetscFunctionBegin; 7788f3b4b4dSDmitry Karpeev if (N < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Total number of subdomains must be 1 or more, got N = %D",N); 7798f3b4b4dSDmitry Karpeev if (pc->setupcalled) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetTotalSubdomains() should be called before calling PCSetUp()."); 7808f3b4b4dSDmitry Karpeev 7812c112581SDmitry Karpeev ierr = PCGASMDestroySubdomains(osm->n,&osm->iis,&osm->ois);CHKERRQ(ierr); 7828f3b4b4dSDmitry Karpeev osm->ois = osm->iis = NULL; 7838f3b4b4dSDmitry Karpeev 7848f3b4b4dSDmitry Karpeev ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr); 7858f3b4b4dSDmitry Karpeev ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);CHKERRQ(ierr); 7868f3b4b4dSDmitry Karpeev osm->N = N; 7878f3b4b4dSDmitry Karpeev osm->n = PETSC_DETERMINE; 7888f3b4b4dSDmitry Karpeev osm->nmax = PETSC_DETERMINE; 7898f3b4b4dSDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 7908f3b4b4dSDmitry Karpeev PetscFunctionReturn(0); 7918f3b4b4dSDmitry Karpeev } 7928f3b4b4dSDmitry Karpeev 7938f3b4b4dSDmitry Karpeev #undef __FUNCT__ 79406b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMSetSubdomains_GASM" 795b541e6a4SDmitry Karpeev static PetscErrorCode PCGASMSetSubdomains_GASM(PC pc,PetscInt n,IS iis[],IS ois[]) 796b1a0a8a3SJed Brown { 797f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 798b1a0a8a3SJed Brown PetscErrorCode ierr; 799b1a0a8a3SJed Brown PetscInt i; 800b1a0a8a3SJed Brown 801b1a0a8a3SJed Brown PetscFunctionBegin; 8028f3b4b4dSDmitry Karpeev if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Each process must have 1 or more subdomains, got n = %D",n); 8038f3b4b4dSDmitry Karpeev if (pc->setupcalled) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetSubdomains() should be called before calling PCSetUp()."); 804b1a0a8a3SJed Brown 8052c112581SDmitry Karpeev ierr = PCGASMDestroySubdomains(osm->n,&osm->iis,&osm->ois);CHKERRQ(ierr); 8068f3b4b4dSDmitry Karpeev osm->iis = osm->ois = NULL; 8078f3b4b4dSDmitry Karpeev osm->n = n; 8088f3b4b4dSDmitry Karpeev osm->N = PETSC_DETERMINE; 8098f3b4b4dSDmitry Karpeev osm->nmax = PETSC_DETERMINE; 810a35b7b57SDmitry Karpeev if (ois) { 811785e854fSJed Brown ierr = PetscMalloc1(n,&osm->ois);CHKERRQ(ierr); 8128f3b4b4dSDmitry Karpeev for (i=0; i<n; i++) { 8138f3b4b4dSDmitry Karpeev ierr = PetscObjectReference((PetscObject)ois[i]);CHKERRQ(ierr); 8148f3b4b4dSDmitry Karpeev osm->ois[i] = ois[i]; 8158f3b4b4dSDmitry Karpeev } 8168f3b4b4dSDmitry Karpeev /* 8178f3b4b4dSDmitry Karpeev Since the user set the outer subdomains, even if nontrivial overlap was requested via PCGASMSetOverlap(), 8188f3b4b4dSDmitry Karpeev it will be ignored. To avoid confusion later on (e.g., when viewing the PC), the overlap size is set to -1. 8198f3b4b4dSDmitry Karpeev */ 820b1a0a8a3SJed Brown osm->overlap = -1; 821a35b7b57SDmitry Karpeev if (!iis) { 822785e854fSJed Brown ierr = PetscMalloc1(n,&osm->iis);CHKERRQ(ierr); 823a35b7b57SDmitry Karpeev for (i=0; i<n; i++) { 8248f3b4b4dSDmitry Karpeev for (i=0; i<n; i++) { 8250e8f3646SDmitry Karpeev ierr = PetscObjectReference((PetscObject)ois[i]);CHKERRQ(ierr); 826a35b7b57SDmitry Karpeev osm->iis[i] = ois[i]; 827a35b7b57SDmitry Karpeev } 828a35b7b57SDmitry Karpeev } 829a35b7b57SDmitry Karpeev } 8308f3b4b4dSDmitry Karpeev } 831a35b7b57SDmitry Karpeev if (iis) { 832785e854fSJed Brown ierr = PetscMalloc1(n,&osm->iis);CHKERRQ(ierr); 8338f3b4b4dSDmitry Karpeev for (i=0; i<n; i++) { 8348f3b4b4dSDmitry Karpeev ierr = PetscObjectReference((PetscObject)iis[i]);CHKERRQ(ierr); 8358f3b4b4dSDmitry Karpeev osm->iis[i] = iis[i]; 8368f3b4b4dSDmitry Karpeev } 837a35b7b57SDmitry Karpeev if (!ois) { 838785e854fSJed Brown ierr = PetscMalloc1(n,&osm->ois);CHKERRQ(ierr); 839a35b7b57SDmitry Karpeev for (i=0; i<n; i++) { 8408f3b4b4dSDmitry Karpeev for (i=0; i<n; i++) { 841a35b7b57SDmitry Karpeev ierr = PetscObjectReference((PetscObject)iis[i]);CHKERRQ(ierr); 842a35b7b57SDmitry Karpeev osm->ois[i] = iis[i]; 843a35b7b57SDmitry Karpeev } 8448f3b4b4dSDmitry Karpeev } 8458f3b4b4dSDmitry Karpeev /* If necessary, osm->ois[i] will be expanded using the requested overlap size in PCSetUp_GASM(). */ 846a35b7b57SDmitry Karpeev } 847a35b7b57SDmitry Karpeev } 8488f3b4b4dSDmitry Karpeev if (iis) osm->user_subdomains = PETSC_TRUE; 849b1a0a8a3SJed Brown PetscFunctionReturn(0); 850b1a0a8a3SJed Brown } 851b1a0a8a3SJed Brown 852b1a0a8a3SJed Brown 853b1a0a8a3SJed Brown #undef __FUNCT__ 854f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetOverlap_GASM" 855f7a08781SBarry Smith static PetscErrorCode PCGASMSetOverlap_GASM(PC pc,PetscInt ovl) 856b1a0a8a3SJed Brown { 857f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 858b1a0a8a3SJed Brown 859b1a0a8a3SJed Brown PetscFunctionBegin; 860ce94432eSBarry Smith if (ovl < 0) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap value requested"); 861ce94432eSBarry Smith if (pc->setupcalled && ovl != osm->overlap) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetOverlap() should be called before PCSetUp()."); 8622fa5cd67SKarl Rupp if (!pc->setupcalled) osm->overlap = ovl; 863b1a0a8a3SJed Brown PetscFunctionReturn(0); 864b1a0a8a3SJed Brown } 865b1a0a8a3SJed Brown 866b1a0a8a3SJed Brown #undef __FUNCT__ 867f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetType_GASM" 868f7a08781SBarry Smith static PetscErrorCode PCGASMSetType_GASM(PC pc,PCGASMType type) 869b1a0a8a3SJed Brown { 870f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 871b1a0a8a3SJed Brown 872b1a0a8a3SJed Brown PetscFunctionBegin; 873b1a0a8a3SJed Brown osm->type = type; 874b1a0a8a3SJed Brown osm->type_set = PETSC_TRUE; 875b1a0a8a3SJed Brown PetscFunctionReturn(0); 876b1a0a8a3SJed Brown } 877b1a0a8a3SJed Brown 878b1a0a8a3SJed Brown #undef __FUNCT__ 879f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetSortIndices_GASM" 880f7a08781SBarry Smith static PetscErrorCode PCGASMSetSortIndices_GASM(PC pc,PetscBool doSort) 881b1a0a8a3SJed Brown { 882f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 883b1a0a8a3SJed Brown 884b1a0a8a3SJed Brown PetscFunctionBegin; 885b1a0a8a3SJed Brown osm->sort_indices = doSort; 886b1a0a8a3SJed Brown PetscFunctionReturn(0); 887b1a0a8a3SJed Brown } 888b1a0a8a3SJed Brown 889b1a0a8a3SJed Brown #undef __FUNCT__ 890f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMGetSubKSP_GASM" 89144fe18e5SDmitry Karpeev /* 8928f3b4b4dSDmitry Karpeev FIXME: This routine might need to be modified now that multiple ranks per subdomain are allowed. 89344fe18e5SDmitry Karpeev In particular, it would upset the global subdomain number calculation. 89444fe18e5SDmitry Karpeev */ 895f7a08781SBarry Smith static PetscErrorCode PCGASMGetSubKSP_GASM(PC pc,PetscInt *n,PetscInt *first,KSP **ksp) 896b1a0a8a3SJed Brown { 897f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 898b1a0a8a3SJed Brown PetscErrorCode ierr; 899b1a0a8a3SJed Brown 900b1a0a8a3SJed Brown PetscFunctionBegin; 901ce94432eSBarry Smith if (osm->n < 1) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Need to call PCSetUP() on PC (or KSPSetUp() on the outer KSP object) before calling here"); 902b1a0a8a3SJed Brown 9032fa5cd67SKarl Rupp if (n) *n = osm->n; 904ab3e923fSDmitry Karpeev if (first) { 905ce94432eSBarry Smith ierr = MPI_Scan(&osm->n,first,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr); 906ab3e923fSDmitry Karpeev *first -= osm->n; 907b1a0a8a3SJed Brown } 908b1a0a8a3SJed Brown if (ksp) { 909b1a0a8a3SJed Brown /* Assume that local solves are now different; not necessarily 91006b43e7eSDmitry Karpeev true, though! This flag is used only for PCView_GASM() */ 911b1a0a8a3SJed Brown *ksp = osm->ksp; 9126a4f0f73SDmitry Karpeev osm->same_subdomain_solvers = PETSC_FALSE; 913b1a0a8a3SJed Brown } 914b1a0a8a3SJed Brown PetscFunctionReturn(0); 915ab3e923fSDmitry Karpeev } /* PCGASMGetSubKSP_GASM() */ 916b1a0a8a3SJed Brown 917b1a0a8a3SJed Brown #undef __FUNCT__ 91806b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMSetSubdomains" 919b1a0a8a3SJed Brown /*@C 92006b43e7eSDmitry Karpeev PCGASMSetSubdomains - Sets the subdomains for this processor 92106b43e7eSDmitry Karpeev for the additive Schwarz preconditioner. 922b1a0a8a3SJed Brown 9238f3b4b4dSDmitry Karpeev Collective on pc 924b1a0a8a3SJed Brown 925b1a0a8a3SJed Brown Input Parameters: 9268f3b4b4dSDmitry Karpeev + pc - the preconditioner object 92706b43e7eSDmitry Karpeev . n - the number of subdomains for this processor 9288f3b4b4dSDmitry Karpeev . iis - the index sets that define the inner subdomains (or NULL for PETSc to determine subdomains) 9298f3b4b4dSDmitry Karpeev - ois - the index sets that define the outer subdomains (or NULL to use the same as iis, or to construct by expanding iis by the requested overlap) 930b1a0a8a3SJed Brown 931b1a0a8a3SJed Brown Notes: 93206b43e7eSDmitry Karpeev The IS indices use the parallel, global numbering of the vector entries. 9336a4f0f73SDmitry Karpeev Inner subdomains are those where the correction is applied. 9346a4f0f73SDmitry Karpeev Outer subdomains are those where the residual necessary to obtain the 9356a4f0f73SDmitry Karpeev corrections is obtained (see PCGASMType for the use of inner/outer subdomains). 9366a4f0f73SDmitry Karpeev Both inner and outer subdomains can extend over several processors. 9376a4f0f73SDmitry Karpeev This processor's portion of a subdomain is known as a local subdomain. 9386a4f0f73SDmitry Karpeev 9396a4f0f73SDmitry Karpeev By default the GASM preconditioner uses 1 (local) subdomain per processor. 9406a4f0f73SDmitry Karpeev 941b1a0a8a3SJed Brown 942b1a0a8a3SJed Brown Level: advanced 943b1a0a8a3SJed Brown 94406b43e7eSDmitry Karpeev .keywords: PC, GASM, set, subdomains, additive Schwarz 945b1a0a8a3SJed Brown 9468f3b4b4dSDmitry Karpeev .seealso: PCGASMSetNumSubdomains(), PCGASMSetOverlap(), PCGASMGetSubKSP(), 94706b43e7eSDmitry Karpeev PCGASMCreateSubdomains2D(), PCGASMGetSubdomains() 948b1a0a8a3SJed Brown @*/ 9496a4f0f73SDmitry Karpeev PetscErrorCode PCGASMSetSubdomains(PC pc,PetscInt n,IS iis[],IS ois[]) 950b1a0a8a3SJed Brown { 9518f3b4b4dSDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 952b1a0a8a3SJed Brown PetscErrorCode ierr; 953b1a0a8a3SJed Brown 954b1a0a8a3SJed Brown PetscFunctionBegin; 955b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 9566a4f0f73SDmitry Karpeev ierr = PetscTryMethod(pc,"PCGASMSetSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,n,iis,ois));CHKERRQ(ierr); 9578f3b4b4dSDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 958b1a0a8a3SJed Brown PetscFunctionReturn(0); 959b1a0a8a3SJed Brown } 960b1a0a8a3SJed Brown 961b1a0a8a3SJed Brown 962b1a0a8a3SJed Brown #undef __FUNCT__ 963f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetOverlap" 964b1a0a8a3SJed Brown /*@ 965f746d493SDmitry Karpeev PCGASMSetOverlap - Sets the overlap between a pair of subdomains for the 966b1a0a8a3SJed Brown additive Schwarz preconditioner. Either all or no processors in the 9678f3b4b4dSDmitry Karpeev pc communicator must call this routine. 968b1a0a8a3SJed Brown 9698f3b4b4dSDmitry Karpeev Logically Collective on pc 970b1a0a8a3SJed Brown 971b1a0a8a3SJed Brown Input Parameters: 972b1a0a8a3SJed Brown + pc - the preconditioner context 9738f3b4b4dSDmitry Karpeev - ovl - the amount of overlap between subdomains (ovl >= 0, default value = 0) 974b1a0a8a3SJed Brown 975b1a0a8a3SJed Brown Options Database Key: 97606b43e7eSDmitry Karpeev . -pc_gasm_overlap <overlap> - Sets overlap 977b1a0a8a3SJed Brown 978b1a0a8a3SJed Brown Notes: 97906b43e7eSDmitry Karpeev By default the GASM preconditioner uses 1 subdomain per processor. To use 9808f3b4b4dSDmitry Karpeev multiple subdomain per perocessor or "straddling" subdomains that intersect 9818f3b4b4dSDmitry Karpeev multiple processors use PCGASMSetSubdomains() (or option -pc_gasm_total_subdomains <n>). 982b1a0a8a3SJed Brown 9838f3b4b4dSDmitry Karpeev The overlap defaults to 0, so if one desires that no additional 984b1a0a8a3SJed Brown overlap be computed beyond what may have been set with a call to 9858f3b4b4dSDmitry Karpeev PCGASMSetSubdomains(), then ovl must be set to be 0. In particular, if one does 9868f3b4b4dSDmitry Karpeev not explicitly set the subdomains in application code, then all overlap would be computed 987f746d493SDmitry Karpeev internally by PETSc, and using an overlap of 0 would result in an GASM 988b1a0a8a3SJed Brown variant that is equivalent to the block Jacobi preconditioner. 989b1a0a8a3SJed Brown 990b1a0a8a3SJed Brown Note that one can define initial index sets with any overlap via 99106b43e7eSDmitry Karpeev PCGASMSetSubdomains(); the routine PCGASMSetOverlap() merely allows 99206b43e7eSDmitry Karpeev PETSc to extend that overlap further, if desired. 993b1a0a8a3SJed Brown 994b1a0a8a3SJed Brown Level: intermediate 995b1a0a8a3SJed Brown 996f746d493SDmitry Karpeev .keywords: PC, GASM, set, overlap 997b1a0a8a3SJed Brown 9988f3b4b4dSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMGetSubKSP(), 99906b43e7eSDmitry Karpeev PCGASMCreateSubdomains2D(), PCGASMGetSubdomains() 1000b1a0a8a3SJed Brown @*/ 10017087cfbeSBarry Smith PetscErrorCode PCGASMSetOverlap(PC pc,PetscInt ovl) 1002b1a0a8a3SJed Brown { 1003b1a0a8a3SJed Brown PetscErrorCode ierr; 10048f3b4b4dSDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 1005b1a0a8a3SJed Brown 1006b1a0a8a3SJed Brown PetscFunctionBegin; 1007b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1008b1a0a8a3SJed Brown PetscValidLogicalCollectiveInt(pc,ovl,2); 1009f746d493SDmitry Karpeev ierr = PetscTryMethod(pc,"PCGASMSetOverlap_C",(PC,PetscInt),(pc,ovl));CHKERRQ(ierr); 10108f3b4b4dSDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 1011b1a0a8a3SJed Brown PetscFunctionReturn(0); 1012b1a0a8a3SJed Brown } 1013b1a0a8a3SJed Brown 1014b1a0a8a3SJed Brown #undef __FUNCT__ 1015f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetType" 1016b1a0a8a3SJed Brown /*@ 1017f746d493SDmitry Karpeev PCGASMSetType - Sets the type of restriction and interpolation used 1018b1a0a8a3SJed Brown for local problems in the additive Schwarz method. 1019b1a0a8a3SJed Brown 1020b1a0a8a3SJed Brown Logically Collective on PC 1021b1a0a8a3SJed Brown 1022b1a0a8a3SJed Brown Input Parameters: 1023b1a0a8a3SJed Brown + pc - the preconditioner context 1024f746d493SDmitry Karpeev - type - variant of GASM, one of 1025b1a0a8a3SJed Brown .vb 1026f746d493SDmitry Karpeev PC_GASM_BASIC - full interpolation and restriction 1027f746d493SDmitry Karpeev PC_GASM_RESTRICT - full restriction, local processor interpolation 1028f746d493SDmitry Karpeev PC_GASM_INTERPOLATE - full interpolation, local processor restriction 1029f746d493SDmitry Karpeev PC_GASM_NONE - local processor restriction and interpolation 1030b1a0a8a3SJed Brown .ve 1031b1a0a8a3SJed Brown 1032b1a0a8a3SJed Brown Options Database Key: 1033f746d493SDmitry Karpeev . -pc_gasm_type [basic,restrict,interpolate,none] - Sets GASM type 1034b1a0a8a3SJed Brown 1035b1a0a8a3SJed Brown Level: intermediate 1036b1a0a8a3SJed Brown 1037f746d493SDmitry Karpeev .keywords: PC, GASM, set, type 1038b1a0a8a3SJed Brown 10398f3b4b4dSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMGetSubKSP(), 1040f746d493SDmitry Karpeev PCGASMCreateSubdomains2D() 1041b1a0a8a3SJed Brown @*/ 10427087cfbeSBarry Smith PetscErrorCode PCGASMSetType(PC pc,PCGASMType type) 1043b1a0a8a3SJed Brown { 1044b1a0a8a3SJed Brown PetscErrorCode ierr; 1045b1a0a8a3SJed Brown 1046b1a0a8a3SJed Brown PetscFunctionBegin; 1047b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1048b1a0a8a3SJed Brown PetscValidLogicalCollectiveEnum(pc,type,2); 1049f746d493SDmitry Karpeev ierr = PetscTryMethod(pc,"PCGASMSetType_C",(PC,PCGASMType),(pc,type));CHKERRQ(ierr); 1050b1a0a8a3SJed Brown PetscFunctionReturn(0); 1051b1a0a8a3SJed Brown } 1052b1a0a8a3SJed Brown 1053b1a0a8a3SJed Brown #undef __FUNCT__ 1054f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetSortIndices" 1055b1a0a8a3SJed Brown /*@ 1056f746d493SDmitry Karpeev PCGASMSetSortIndices - Determines whether subdomain indices are sorted. 1057b1a0a8a3SJed Brown 1058b1a0a8a3SJed Brown Logically Collective on PC 1059b1a0a8a3SJed Brown 1060b1a0a8a3SJed Brown Input Parameters: 1061b1a0a8a3SJed Brown + pc - the preconditioner context 1062b1a0a8a3SJed Brown - doSort - sort the subdomain indices 1063b1a0a8a3SJed Brown 1064b1a0a8a3SJed Brown Level: intermediate 1065b1a0a8a3SJed Brown 1066f746d493SDmitry Karpeev .keywords: PC, GASM, set, type 1067b1a0a8a3SJed Brown 10688f3b4b4dSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMGetSubKSP(), 1069f746d493SDmitry Karpeev PCGASMCreateSubdomains2D() 1070b1a0a8a3SJed Brown @*/ 10717087cfbeSBarry Smith PetscErrorCode PCGASMSetSortIndices(PC pc,PetscBool doSort) 1072b1a0a8a3SJed Brown { 1073b1a0a8a3SJed Brown PetscErrorCode ierr; 1074b1a0a8a3SJed Brown 1075b1a0a8a3SJed Brown PetscFunctionBegin; 1076b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1077b1a0a8a3SJed Brown PetscValidLogicalCollectiveBool(pc,doSort,2); 1078f746d493SDmitry Karpeev ierr = PetscTryMethod(pc,"PCGASMSetSortIndices_C",(PC,PetscBool),(pc,doSort));CHKERRQ(ierr); 1079b1a0a8a3SJed Brown PetscFunctionReturn(0); 1080b1a0a8a3SJed Brown } 1081b1a0a8a3SJed Brown 1082b1a0a8a3SJed Brown #undef __FUNCT__ 1083f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMGetSubKSP" 1084b1a0a8a3SJed Brown /*@C 1085f746d493SDmitry Karpeev PCGASMGetSubKSP - Gets the local KSP contexts for all blocks on 1086b1a0a8a3SJed Brown this processor. 1087b1a0a8a3SJed Brown 1088b1a0a8a3SJed Brown Collective on PC iff first_local is requested 1089b1a0a8a3SJed Brown 1090b1a0a8a3SJed Brown Input Parameter: 1091b1a0a8a3SJed Brown . pc - the preconditioner context 1092b1a0a8a3SJed Brown 1093b1a0a8a3SJed Brown Output Parameters: 10940298fd71SBarry Smith + n_local - the number of blocks on this processor or NULL 10950298fd71SBarry Smith . first_local - the global number of the first block on this processor or NULL, 10960298fd71SBarry Smith all processors must request or all must pass NULL 1097b1a0a8a3SJed Brown - ksp - the array of KSP contexts 1098b1a0a8a3SJed Brown 1099b1a0a8a3SJed Brown Note: 1100f746d493SDmitry Karpeev After PCGASMGetSubKSP() the array of KSPes is not to be freed 1101b1a0a8a3SJed Brown 1102b1a0a8a3SJed Brown Currently for some matrix implementations only 1 block per processor 1103b1a0a8a3SJed Brown is supported. 1104b1a0a8a3SJed Brown 1105f746d493SDmitry Karpeev You must call KSPSetUp() before calling PCGASMGetSubKSP(). 1106b1a0a8a3SJed Brown 1107b1a0a8a3SJed Brown Level: advanced 1108b1a0a8a3SJed Brown 1109f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, get, sub, KSP, context 1110b1a0a8a3SJed Brown 11118f3b4b4dSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMSetOverlap(), 1112f746d493SDmitry Karpeev PCGASMCreateSubdomains2D(), 1113b1a0a8a3SJed Brown @*/ 11147087cfbeSBarry Smith PetscErrorCode PCGASMGetSubKSP(PC pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp[]) 1115b1a0a8a3SJed Brown { 1116b1a0a8a3SJed Brown PetscErrorCode ierr; 1117b1a0a8a3SJed Brown 1118b1a0a8a3SJed Brown PetscFunctionBegin; 1119b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1120f746d493SDmitry Karpeev ierr = PetscUseMethod(pc,"PCGASMGetSubKSP_C",(PC,PetscInt*,PetscInt*,KSP **),(pc,n_local,first_local,ksp));CHKERRQ(ierr); 1121b1a0a8a3SJed Brown PetscFunctionReturn(0); 1122b1a0a8a3SJed Brown } 1123b1a0a8a3SJed Brown 1124b1a0a8a3SJed Brown /* -------------------------------------------------------------------------------------*/ 1125b1a0a8a3SJed Brown /*MC 1126f746d493SDmitry Karpeev PCGASM - Use the (restricted) additive Schwarz method, each block is (approximately) solved with 1127b1a0a8a3SJed Brown its own KSP object. 1128b1a0a8a3SJed Brown 1129b1a0a8a3SJed Brown Options Database Keys: 11308f3b4b4dSDmitry Karpeev + -pc_gasm_total_subdomains <n> - Sets total number of local subdomains to be distributed among processors 113106b43e7eSDmitry Karpeev . -pc_gasm_view_subdomains - activates the printing of subdomain indices in PCView(), -ksp_view or -snes_view 113206b43e7eSDmitry Karpeev . -pc_gasm_print_subdomains - activates the printing of subdomain indices in PCSetUp() 113306b43e7eSDmitry Karpeev . -pc_gasm_overlap <ovl> - Sets overlap by which to (automatically) extend local subdomains 1134f746d493SDmitry Karpeev - -pc_gasm_type [basic,restrict,interpolate,none] - Sets GASM type 1135b1a0a8a3SJed Brown 1136b1a0a8a3SJed Brown IMPORTANT: If you run with, for example, 3 blocks on 1 processor or 3 blocks on 3 processors you 1137f746d493SDmitry Karpeev will get a different convergence rate due to the default option of -pc_gasm_type restrict. Use 1138f746d493SDmitry Karpeev -pc_gasm_type basic to use the standard GASM. 1139b1a0a8a3SJed Brown 1140b1a0a8a3SJed Brown Notes: Each processor can have one or more blocks, but a block cannot be shared by more 1141b1a0a8a3SJed Brown than one processor. Defaults to one block per processor. 1142b1a0a8a3SJed Brown 1143b1a0a8a3SJed Brown To set options on the solvers for each block append -sub_ to all the KSP, and PC 1144b1a0a8a3SJed Brown options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly 1145b1a0a8a3SJed Brown 1146f746d493SDmitry Karpeev To set the options on the solvers separate for each block call PCGASMGetSubKSP() 1147b1a0a8a3SJed Brown and set the options directly on the resulting KSP object (you can access its PC 1148b1a0a8a3SJed Brown with KSPGetPC()) 1149b1a0a8a3SJed Brown 1150b1a0a8a3SJed Brown 1151b1a0a8a3SJed Brown Level: beginner 1152b1a0a8a3SJed Brown 1153b1a0a8a3SJed Brown Concepts: additive Schwarz method 1154b1a0a8a3SJed Brown 1155b1a0a8a3SJed Brown References: 1156b1a0a8a3SJed Brown An additive variant of the Schwarz alternating method for the case of many subregions 1157b1a0a8a3SJed Brown M Dryja, OB Widlund - Courant Institute, New York University Technical report 1158b1a0a8a3SJed Brown 1159b1a0a8a3SJed Brown Domain Decompositions: Parallel Multilevel Methods for Elliptic Partial Differential Equations, 1160b1a0a8a3SJed Brown Barry Smith, Petter Bjorstad, and William Gropp, Cambridge University Press, ISBN 0-521-49589-X. 1161b1a0a8a3SJed Brown 1162b1a0a8a3SJed Brown .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 116349517cdeSBarry Smith PCBJACOBI, PCGASMGetSubKSP(), PCGASMSetSubdomains(), 11648f3b4b4dSDmitry Karpeev PCSetModifySubmatrices(), PCGASMSetOverlap(), PCGASMSetType() 1165b1a0a8a3SJed Brown 1166b1a0a8a3SJed Brown M*/ 1167b1a0a8a3SJed Brown 1168b1a0a8a3SJed Brown #undef __FUNCT__ 1169f746d493SDmitry Karpeev #define __FUNCT__ "PCCreate_GASM" 11708cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_GASM(PC pc) 1171b1a0a8a3SJed Brown { 1172b1a0a8a3SJed Brown PetscErrorCode ierr; 1173f746d493SDmitry Karpeev PC_GASM *osm; 1174b1a0a8a3SJed Brown 1175b1a0a8a3SJed Brown PetscFunctionBegin; 1176b00a9115SJed Brown ierr = PetscNewLog(pc,&osm);CHKERRQ(ierr); 11772fa5cd67SKarl Rupp 11788f3b4b4dSDmitry Karpeev osm->N = PETSC_DETERMINE; 117906b43e7eSDmitry Karpeev osm->n = PETSC_DECIDE; 11808f3b4b4dSDmitry Karpeev osm->nmax = PETSC_DETERMINE; 11818f3b4b4dSDmitry Karpeev osm->overlap = 0; 1182b1a0a8a3SJed Brown osm->ksp = 0; 11836a4f0f73SDmitry Karpeev osm->gorestriction = 0; 11846a4f0f73SDmitry Karpeev osm->girestriction = 0; 1185ab3e923fSDmitry Karpeev osm->gx = 0; 1186ab3e923fSDmitry Karpeev osm->gy = 0; 1187b1a0a8a3SJed Brown osm->x = 0; 1188b1a0a8a3SJed Brown osm->y = 0; 11896a4f0f73SDmitry Karpeev osm->ois = 0; 11906a4f0f73SDmitry Karpeev osm->iis = 0; 1191b1a0a8a3SJed Brown osm->pmat = 0; 1192f746d493SDmitry Karpeev osm->type = PC_GASM_RESTRICT; 11936a4f0f73SDmitry Karpeev osm->same_subdomain_solvers = PETSC_TRUE; 1194b1a0a8a3SJed Brown osm->sort_indices = PETSC_TRUE; 1195d709ea83SDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 1196b1a0a8a3SJed Brown 1197b1a0a8a3SJed Brown pc->data = (void*)osm; 1198f746d493SDmitry Karpeev pc->ops->apply = PCApply_GASM; 1199f746d493SDmitry Karpeev pc->ops->applytranspose = PCApplyTranspose_GASM; 1200f746d493SDmitry Karpeev pc->ops->setup = PCSetUp_GASM; 1201a06653b4SBarry Smith pc->ops->reset = PCReset_GASM; 1202f746d493SDmitry Karpeev pc->ops->destroy = PCDestroy_GASM; 1203f746d493SDmitry Karpeev pc->ops->setfromoptions = PCSetFromOptions_GASM; 1204f746d493SDmitry Karpeev pc->ops->setuponblocks = PCSetUpOnBlocks_GASM; 1205f746d493SDmitry Karpeev pc->ops->view = PCView_GASM; 1206b1a0a8a3SJed Brown pc->ops->applyrichardson = 0; 1207b1a0a8a3SJed Brown 1208bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetSubdomains_C",PCGASMSetSubdomains_GASM);CHKERRQ(ierr); 1209bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetOverlap_C",PCGASMSetOverlap_GASM);CHKERRQ(ierr); 1210bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetType_C",PCGASMSetType_GASM);CHKERRQ(ierr); 1211bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetSortIndices_C",PCGASMSetSortIndices_GASM);CHKERRQ(ierr); 1212bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGASMGetSubKSP_C",PCGASMGetSubKSP_GASM);CHKERRQ(ierr); 1213b1a0a8a3SJed Brown PetscFunctionReturn(0); 1214b1a0a8a3SJed Brown } 1215b1a0a8a3SJed Brown 1216b1a0a8a3SJed Brown 1217b1a0a8a3SJed Brown #undef __FUNCT__ 121806b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMCreateLocalSubdomains" 12198f3b4b4dSDmitry Karpeev PetscErrorCode PCGASMCreateLocalSubdomains(Mat A, PetscInt nloc, IS *iis[]) 1220b1a0a8a3SJed Brown { 1221b1a0a8a3SJed Brown MatPartitioning mpart; 1222b1a0a8a3SJed Brown const char *prefix; 122311bd1e4dSLisandro Dalcin PetscErrorCode (*f)(Mat,MatReuse,Mat*); 1224b1a0a8a3SJed Brown PetscMPIInt size; 1225b1a0a8a3SJed Brown PetscInt i,j,rstart,rend,bs; 122611bd1e4dSLisandro Dalcin PetscBool isbaij = PETSC_FALSE,foundpart = PETSC_FALSE; 12270298fd71SBarry Smith Mat Ad = NULL, adj; 1228b1a0a8a3SJed Brown IS ispart,isnumb,*is; 1229b1a0a8a3SJed Brown PetscErrorCode ierr; 1230b1a0a8a3SJed Brown 1231b1a0a8a3SJed Brown PetscFunctionBegin; 12328f3b4b4dSDmitry Karpeev if (nloc < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of local subdomains must > 0, got nloc = %D",nloc); 1233b1a0a8a3SJed Brown 1234b1a0a8a3SJed Brown /* Get prefix, row distribution, and block size */ 1235b1a0a8a3SJed Brown ierr = MatGetOptionsPrefix(A,&prefix);CHKERRQ(ierr); 1236b1a0a8a3SJed Brown ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 1237b1a0a8a3SJed Brown ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr); 1238b1a0a8a3SJed Brown 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); 1239b1a0a8a3SJed Brown 1240b1a0a8a3SJed Brown /* Get diagonal block from matrix if possible */ 1241ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr); 12420005d66cSJed Brown ierr = PetscObjectQueryFunction((PetscObject)A,"MatGetDiagonalBlock_C",&f);CHKERRQ(ierr); 1243b1a0a8a3SJed Brown if (f) { 124411bd1e4dSLisandro Dalcin ierr = MatGetDiagonalBlock(A,&Ad);CHKERRQ(ierr); 1245b1a0a8a3SJed Brown } else if (size == 1) { 124611bd1e4dSLisandro Dalcin Ad = A; 1247b1a0a8a3SJed Brown } 1248b1a0a8a3SJed Brown if (Ad) { 1249251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)Ad,MATSEQBAIJ,&isbaij);CHKERRQ(ierr); 1250251f4c67SDmitry Karpeev if (!isbaij) {ierr = PetscObjectTypeCompare((PetscObject)Ad,MATSEQSBAIJ,&isbaij);CHKERRQ(ierr);} 1251b1a0a8a3SJed Brown } 12528f3b4b4dSDmitry Karpeev if (Ad && nloc > 1) { 1253b1a0a8a3SJed Brown PetscBool match,done; 1254b1a0a8a3SJed Brown /* Try to setup a good matrix partitioning if available */ 1255b1a0a8a3SJed Brown ierr = MatPartitioningCreate(PETSC_COMM_SELF,&mpart);CHKERRQ(ierr); 1256b1a0a8a3SJed Brown ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr); 1257b1a0a8a3SJed Brown ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr); 1258251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGCURRENT,&match);CHKERRQ(ierr); 1259b1a0a8a3SJed Brown if (!match) { 1260251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGSQUARE,&match);CHKERRQ(ierr); 1261b1a0a8a3SJed Brown } 1262b1a0a8a3SJed Brown if (!match) { /* assume a "good" partitioner is available */ 12631a83f524SJed Brown PetscInt na; 12641a83f524SJed Brown const PetscInt *ia,*ja; 1265b1a0a8a3SJed Brown ierr = MatGetRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr); 1266b1a0a8a3SJed Brown if (done) { 1267b1a0a8a3SJed Brown /* Build adjacency matrix by hand. Unfortunately a call to 1268b1a0a8a3SJed Brown MatConvert(Ad,MATMPIADJ,MAT_INITIAL_MATRIX,&adj) will 1269b1a0a8a3SJed Brown remove the block-aij structure and we cannot expect 1270b1a0a8a3SJed Brown MatPartitioning to split vertices as we need */ 12711a83f524SJed Brown PetscInt i,j,len,nnz,cnt,*iia=0,*jja=0; 12721a83f524SJed Brown const PetscInt *row; 1273b1a0a8a3SJed Brown nnz = 0; 1274b1a0a8a3SJed Brown for (i=0; i<na; i++) { /* count number of nonzeros */ 1275b1a0a8a3SJed Brown len = ia[i+1] - ia[i]; 1276b1a0a8a3SJed Brown row = ja + ia[i]; 1277b1a0a8a3SJed Brown for (j=0; j<len; j++) { 1278b1a0a8a3SJed Brown if (row[j] == i) { /* don't count diagonal */ 1279b1a0a8a3SJed Brown len--; break; 1280b1a0a8a3SJed Brown } 1281b1a0a8a3SJed Brown } 1282b1a0a8a3SJed Brown nnz += len; 1283b1a0a8a3SJed Brown } 1284854ce69bSBarry Smith ierr = PetscMalloc1(na+1,&iia);CHKERRQ(ierr); 1285854ce69bSBarry Smith ierr = PetscMalloc1(nnz,&jja);CHKERRQ(ierr); 1286b1a0a8a3SJed Brown nnz = 0; 1287b1a0a8a3SJed Brown iia[0] = 0; 1288b1a0a8a3SJed Brown for (i=0; i<na; i++) { /* fill adjacency */ 1289b1a0a8a3SJed Brown cnt = 0; 1290b1a0a8a3SJed Brown len = ia[i+1] - ia[i]; 1291b1a0a8a3SJed Brown row = ja + ia[i]; 1292b1a0a8a3SJed Brown for (j=0; j<len; j++) { 12932fa5cd67SKarl Rupp if (row[j] != i) jja[nnz+cnt++] = row[j]; /* if not diagonal */ 1294b1a0a8a3SJed Brown } 1295b1a0a8a3SJed Brown nnz += cnt; 1296b1a0a8a3SJed Brown iia[i+1] = nnz; 1297b1a0a8a3SJed Brown } 1298b1a0a8a3SJed Brown /* Partitioning of the adjacency matrix */ 12990298fd71SBarry Smith ierr = MatCreateMPIAdj(PETSC_COMM_SELF,na,na,iia,jja,NULL,&adj);CHKERRQ(ierr); 1300b1a0a8a3SJed Brown ierr = MatPartitioningSetAdjacency(mpart,adj);CHKERRQ(ierr); 13018f3b4b4dSDmitry Karpeev ierr = MatPartitioningSetNParts(mpart,nloc);CHKERRQ(ierr); 1302b1a0a8a3SJed Brown ierr = MatPartitioningApply(mpart,&ispart);CHKERRQ(ierr); 1303b1a0a8a3SJed Brown ierr = ISPartitioningToNumbering(ispart,&isnumb);CHKERRQ(ierr); 13046bf464f9SBarry Smith ierr = MatDestroy(&adj);CHKERRQ(ierr); 1305b1a0a8a3SJed Brown foundpart = PETSC_TRUE; 1306b1a0a8a3SJed Brown } 1307b1a0a8a3SJed Brown ierr = MatRestoreRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr); 1308b1a0a8a3SJed Brown } 13096bf464f9SBarry Smith ierr = MatPartitioningDestroy(&mpart);CHKERRQ(ierr); 1310b1a0a8a3SJed Brown } 13118f3b4b4dSDmitry Karpeev ierr = PetscMalloc1(nloc,&is);CHKERRQ(ierr); 1312b1a0a8a3SJed Brown if (!foundpart) { 1313b1a0a8a3SJed Brown 1314b1a0a8a3SJed Brown /* Partitioning by contiguous chunks of rows */ 1315b1a0a8a3SJed Brown 1316b1a0a8a3SJed Brown PetscInt mbs = (rend-rstart)/bs; 1317b1a0a8a3SJed Brown PetscInt start = rstart; 13188f3b4b4dSDmitry Karpeev for (i=0; i<nloc; i++) { 13198f3b4b4dSDmitry Karpeev PetscInt count = (mbs/nloc + ((mbs % nloc) > i)) * bs; 1320b1a0a8a3SJed Brown ierr = ISCreateStride(PETSC_COMM_SELF,count,start,1,&is[i]);CHKERRQ(ierr); 1321b1a0a8a3SJed Brown start += count; 1322b1a0a8a3SJed Brown } 1323b1a0a8a3SJed Brown 1324b1a0a8a3SJed Brown } else { 1325b1a0a8a3SJed Brown 1326b1a0a8a3SJed Brown /* Partitioning by adjacency of diagonal block */ 1327b1a0a8a3SJed Brown 1328b1a0a8a3SJed Brown const PetscInt *numbering; 1329b1a0a8a3SJed Brown PetscInt *count,nidx,*indices,*newidx,start=0; 1330b1a0a8a3SJed Brown /* Get node count in each partition */ 13318f3b4b4dSDmitry Karpeev ierr = PetscMalloc1(nloc,&count);CHKERRQ(ierr); 13328f3b4b4dSDmitry Karpeev ierr = ISPartitioningCount(ispart,nloc,count);CHKERRQ(ierr); 1333b1a0a8a3SJed Brown if (isbaij && bs > 1) { /* adjust for the block-aij case */ 13348f3b4b4dSDmitry Karpeev for (i=0; i<nloc; i++) count[i] *= bs; 1335b1a0a8a3SJed Brown } 1336b1a0a8a3SJed Brown /* Build indices from node numbering */ 1337b1a0a8a3SJed Brown ierr = ISGetLocalSize(isnumb,&nidx);CHKERRQ(ierr); 1338785e854fSJed Brown ierr = PetscMalloc1(nidx,&indices);CHKERRQ(ierr); 1339b1a0a8a3SJed Brown for (i=0; i<nidx; i++) indices[i] = i; /* needs to be initialized */ 1340b1a0a8a3SJed Brown ierr = ISGetIndices(isnumb,&numbering);CHKERRQ(ierr); 1341b1a0a8a3SJed Brown ierr = PetscSortIntWithPermutation(nidx,numbering,indices);CHKERRQ(ierr); 1342b1a0a8a3SJed Brown ierr = ISRestoreIndices(isnumb,&numbering);CHKERRQ(ierr); 1343b1a0a8a3SJed Brown if (isbaij && bs > 1) { /* adjust for the block-aij case */ 1344785e854fSJed Brown ierr = PetscMalloc1(nidx*bs,&newidx);CHKERRQ(ierr); 13452fa5cd67SKarl Rupp for (i=0; i<nidx; i++) { 13462fa5cd67SKarl Rupp for (j=0; j<bs; j++) newidx[i*bs+j] = indices[i]*bs + j; 13472fa5cd67SKarl Rupp } 1348b1a0a8a3SJed Brown ierr = PetscFree(indices);CHKERRQ(ierr); 1349b1a0a8a3SJed Brown nidx *= bs; 1350b1a0a8a3SJed Brown indices = newidx; 1351b1a0a8a3SJed Brown } 1352b1a0a8a3SJed Brown /* Shift to get global indices */ 1353b1a0a8a3SJed Brown for (i=0; i<nidx; i++) indices[i] += rstart; 1354b1a0a8a3SJed Brown 1355b1a0a8a3SJed Brown /* Build the index sets for each block */ 13568f3b4b4dSDmitry Karpeev for (i=0; i<nloc; i++) { 1357b1a0a8a3SJed Brown ierr = ISCreateGeneral(PETSC_COMM_SELF,count[i],&indices[start],PETSC_COPY_VALUES,&is[i]);CHKERRQ(ierr); 1358b1a0a8a3SJed Brown ierr = ISSort(is[i]);CHKERRQ(ierr); 1359b1a0a8a3SJed Brown start += count[i]; 1360b1a0a8a3SJed Brown } 1361b1a0a8a3SJed Brown 1362302440fdSBarry Smith ierr = PetscFree(count);CHKERRQ(ierr); 1363302440fdSBarry Smith ierr = PetscFree(indices);CHKERRQ(ierr); 1364fcfd50ebSBarry Smith ierr = ISDestroy(&isnumb);CHKERRQ(ierr); 1365fcfd50ebSBarry Smith ierr = ISDestroy(&ispart);CHKERRQ(ierr); 1366b1a0a8a3SJed Brown } 13676a4f0f73SDmitry Karpeev *iis = is; 13688f3b4b4dSDmitry Karpeev PetscFunctionReturn(0); 13698f3b4b4dSDmitry Karpeev } 13708f3b4b4dSDmitry Karpeev 13718f3b4b4dSDmitry Karpeev #undef __FUNCT__ 13728f3b4b4dSDmitry Karpeev #define __FUNCT__ "PCGASMCreateStraddlingSubdomains" 1373b541e6a4SDmitry Karpeev PETSC_INTERN PetscErrorCode PCGASMCreateStraddlingSubdomains(Mat A,PetscInt N,PetscInt *n,IS *iis[]) 13748f3b4b4dSDmitry Karpeev { 13758f3b4b4dSDmitry Karpeev PetscErrorCode ierr; 13768f3b4b4dSDmitry Karpeev 13778f3b4b4dSDmitry Karpeev PetscFunctionBegin; 13788f3b4b4dSDmitry Karpeev ierr = MatSubdomainsCreateCoalesce(A,N,n,iis);CHKERRQ(ierr); 13798f3b4b4dSDmitry Karpeev PetscFunctionReturn(0); 13808f3b4b4dSDmitry Karpeev } 13818f3b4b4dSDmitry Karpeev 13828f3b4b4dSDmitry Karpeev 13838f3b4b4dSDmitry Karpeev 13848f3b4b4dSDmitry Karpeev #undef __FUNCT__ 13858f3b4b4dSDmitry Karpeev #define __FUNCT__ "PCGASMCreateSubdomains" 13868f3b4b4dSDmitry Karpeev /*@C 13878f3b4b4dSDmitry Karpeev PCGASMCreateSubdomains - Creates n index sets defining n nonoverlapping subdomains for the additive 13888f3b4b4dSDmitry Karpeev Schwarz preconditioner for a any problem based on its matrix. 13898f3b4b4dSDmitry Karpeev 13908f3b4b4dSDmitry Karpeev Collective 13918f3b4b4dSDmitry Karpeev 13928f3b4b4dSDmitry Karpeev Input Parameters: 13938f3b4b4dSDmitry Karpeev + A - The global matrix operator 13948f3b4b4dSDmitry Karpeev - N - the number of global subdomains requested 13958f3b4b4dSDmitry Karpeev 13968f3b4b4dSDmitry Karpeev Output Parameters: 13978f3b4b4dSDmitry Karpeev + n - the number of subdomains created on this processor 13988f3b4b4dSDmitry Karpeev - iis - the array of index sets defining the local inner subdomains (on which the correction is applied) 13998f3b4b4dSDmitry Karpeev 14008f3b4b4dSDmitry Karpeev Level: advanced 14018f3b4b4dSDmitry Karpeev 14028f3b4b4dSDmitry Karpeev Note: When N >= A's communicator size, each subdomain is local -- contained within a single processor. 14038f3b4b4dSDmitry Karpeev When N < size, the subdomains are 'straddling' (processor boundaries) and are no longer local. 14048f3b4b4dSDmitry Karpeev The resulting subdomains can be use in PCGASMSetSubdomains(pc,n,iss,NULL). The overlapping 14058f3b4b4dSDmitry Karpeev outer subdomains will be automatically generated from these according to the requested amount of 14068f3b4b4dSDmitry Karpeev overlap; this is currently supported only with local subdomains. 14078f3b4b4dSDmitry Karpeev 14088f3b4b4dSDmitry Karpeev 14098f3b4b4dSDmitry Karpeev .keywords: PC, GASM, additive Schwarz, create, subdomains, unstructured grid 14108f3b4b4dSDmitry Karpeev 14118f3b4b4dSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMDestroySubdomains() 14128f3b4b4dSDmitry Karpeev @*/ 1413b541e6a4SDmitry Karpeev PetscErrorCode PCGASMCreateSubdomains(Mat A,PetscInt N,PetscInt *n,IS *iis[]) 14148f3b4b4dSDmitry Karpeev { 14158f3b4b4dSDmitry Karpeev PetscMPIInt size; 14168f3b4b4dSDmitry Karpeev PetscErrorCode ierr; 14178f3b4b4dSDmitry Karpeev 14188f3b4b4dSDmitry Karpeev PetscFunctionBegin; 14198f3b4b4dSDmitry Karpeev PetscValidHeaderSpecific(A,MAT_CLASSID,1); 14208f3b4b4dSDmitry Karpeev PetscValidPointer(iis,4); 14218f3b4b4dSDmitry Karpeev 14228f3b4b4dSDmitry Karpeev if (N < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Number of subdomains must be > 0, N = %D",N); 14238f3b4b4dSDmitry Karpeev ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr); 14248f3b4b4dSDmitry Karpeev if (N >= size) { 14258f3b4b4dSDmitry Karpeev *n = N/size + (N%size); 14268f3b4b4dSDmitry Karpeev ierr = PCGASMCreateLocalSubdomains(A,*n,iis);CHKERRQ(ierr); 14276a4f0f73SDmitry Karpeev } else { 14288f3b4b4dSDmitry Karpeev ierr = PCGASMCreateStraddlingSubdomains(A,N,n,iis);CHKERRQ(ierr); 14296a4f0f73SDmitry Karpeev } 1430b1a0a8a3SJed Brown PetscFunctionReturn(0); 1431b1a0a8a3SJed Brown } 1432b1a0a8a3SJed Brown 1433b1a0a8a3SJed Brown #undef __FUNCT__ 1434f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMDestroySubdomains" 1435b1a0a8a3SJed Brown /*@C 1436f746d493SDmitry Karpeev PCGASMDestroySubdomains - Destroys the index sets created with 14378f3b4b4dSDmitry Karpeev PCGASMCreateSubdomains() or PCGASMCreateSubdomains2D. Should be 143806b43e7eSDmitry Karpeev called after setting subdomains with PCGASMSetSubdomains(). 1439b1a0a8a3SJed Brown 1440b1a0a8a3SJed Brown Collective 1441b1a0a8a3SJed Brown 1442b1a0a8a3SJed Brown Input Parameters: 1443b1a0a8a3SJed Brown + n - the number of index sets 14446a4f0f73SDmitry Karpeev . iis - the array of inner subdomains, 14450298fd71SBarry Smith - ois - the array of outer subdomains, can be NULL 1446b1a0a8a3SJed Brown 14476a4f0f73SDmitry Karpeev Level: intermediate 14486a4f0f73SDmitry Karpeev 14496a4f0f73SDmitry Karpeev Notes: this is merely a convenience subroutine that walks each list, 14502c112581SDmitry Karpeev destroys each IS on the list, and then frees the list. At the end the 14512c112581SDmitry Karpeev list pointers are set to NULL. 1452b1a0a8a3SJed Brown 1453f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, create, subdomains, unstructured grid 1454b1a0a8a3SJed Brown 14558f3b4b4dSDmitry Karpeev .seealso: PCGASMCreateSubdomains(), PCGASMSetSubdomains() 1456b1a0a8a3SJed Brown @*/ 14572c112581SDmitry Karpeev PetscErrorCode PCGASMDestroySubdomains(PetscInt n,IS **iis,IS **ois) 1458b1a0a8a3SJed Brown { 1459b1a0a8a3SJed Brown PetscInt i; 1460b1a0a8a3SJed Brown PetscErrorCode ierr; 14615fd66863SKarl Rupp 1462b1a0a8a3SJed Brown PetscFunctionBegin; 1463a35b7b57SDmitry Karpeev if (n <= 0) PetscFunctionReturn(0); 14646a4f0f73SDmitry Karpeev if (ois) { 14652c112581SDmitry Karpeev PetscValidPointer(ois,3); 14662c112581SDmitry Karpeev if (*ois) { 14672c112581SDmitry Karpeev PetscValidPointer(*ois,3); 1468a35b7b57SDmitry Karpeev for (i=0; i<n; i++) { 14692c112581SDmitry Karpeev ierr = ISDestroy(&(*ois)[i]);CHKERRQ(ierr); 1470a35b7b57SDmitry Karpeev } 1471135757f6SDmitry Karpeev ierr = PetscFree((*ois));CHKERRQ(ierr); 14722c112581SDmitry Karpeev } 1473b1a0a8a3SJed Brown } 1474b9d0fdaaSFande Kong if (iis) { 1475b9d0fdaaSFande Kong PetscValidPointer(iis,2); 1476b9d0fdaaSFande Kong if (*iis) { 1477b9d0fdaaSFande Kong PetscValidPointer(*iis,2); 1478b9d0fdaaSFande Kong for (i=0; i<n; i++) { 1479b9d0fdaaSFande Kong ierr = ISDestroy(&(*iis)[i]);CHKERRQ(ierr); 1480b9d0fdaaSFande Kong } 1481b9d0fdaaSFande Kong ierr = PetscFree((*iis));CHKERRQ(ierr); 1482b9d0fdaaSFande Kong } 1483b9d0fdaaSFande Kong } 1484b1a0a8a3SJed Brown PetscFunctionReturn(0); 1485b1a0a8a3SJed Brown } 1486b1a0a8a3SJed Brown 1487af538404SDmitry Karpeev 1488af538404SDmitry Karpeev #define PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,xleft_loc,ylow_loc,xright_loc,yhigh_loc,n) \ 1489af538404SDmitry Karpeev { \ 1490af538404SDmitry Karpeev PetscInt first_row = first/M, last_row = last/M+1; \ 1491af538404SDmitry Karpeev /* \ 1492af538404SDmitry Karpeev Compute ylow_loc and yhigh_loc so that (ylow_loc,xleft) and (yhigh_loc,xright) are the corners \ 1493af538404SDmitry Karpeev of the bounding box of the intersection of the subdomain with the local ownership range (local \ 1494af538404SDmitry Karpeev subdomain). \ 1495af538404SDmitry Karpeev Also compute xleft_loc and xright_loc as the lower and upper bounds on the first and last rows \ 1496af538404SDmitry Karpeev of the intersection. \ 1497af538404SDmitry Karpeev */ \ 1498af538404SDmitry Karpeev /* ylow_loc is the grid row containing the first element of the local sumbdomain */ \ 1499eec7e3faSDmitry Karpeev *ylow_loc = PetscMax(first_row,ylow); \ 1500af538404SDmitry Karpeev /* xleft_loc is the offset of first element of the local subdomain within its grid row (might actually be outside the local subdomain) */ \ 1501eec7e3faSDmitry Karpeev *xleft_loc = *ylow_loc==first_row ? PetscMax(first%M,xleft) : xleft; \ 1502af538404SDmitry Karpeev /* yhigh_loc is the grid row above the last local subdomain element */ \ 1503eec7e3faSDmitry Karpeev *yhigh_loc = PetscMin(last_row,yhigh); \ 1504af538404SDmitry Karpeev /* xright is the offset of the end of the local subdomain within its grid row (might actually be outside the local subdomain) */ \ 1505eec7e3faSDmitry Karpeev *xright_loc = *yhigh_loc==last_row ? PetscMin(xright,last%M) : xright; \ 1506af538404SDmitry Karpeev /* Now compute the size of the local subdomain n. */ \ 1507c3518bceSDmitry Karpeev *n = 0; \ 1508eec7e3faSDmitry Karpeev if (*ylow_loc < *yhigh_loc) { \ 1509af538404SDmitry Karpeev PetscInt width = xright-xleft; \ 1510eec7e3faSDmitry Karpeev *n += width*(*yhigh_loc-*ylow_loc-1); \ 1511eec7e3faSDmitry Karpeev *n += PetscMin(PetscMax(*xright_loc-xleft,0),width); \ 1512eec7e3faSDmitry Karpeev *n -= PetscMin(PetscMax(*xleft_loc-xleft,0), width); \ 1513af538404SDmitry Karpeev } \ 1514af538404SDmitry Karpeev } 1515af538404SDmitry Karpeev 1516c3518bceSDmitry Karpeev 1517eec7e3faSDmitry Karpeev 1518eec7e3faSDmitry Karpeev #undef __FUNCT__ 1519f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMCreateSubdomains2D" 1520b1a0a8a3SJed Brown /*@ 1521f746d493SDmitry Karpeev PCGASMCreateSubdomains2D - Creates the index sets for the overlapping Schwarz 1522b1a0a8a3SJed Brown preconditioner for a two-dimensional problem on a regular grid. 1523b1a0a8a3SJed Brown 1524af538404SDmitry Karpeev Collective 1525b1a0a8a3SJed Brown 1526b1a0a8a3SJed Brown Input Parameters: 152706b43e7eSDmitry Karpeev + M, N - the global number of grid points in the x and y directions 1528af538404SDmitry Karpeev . Mdomains, Ndomains - the global number of subdomains in the x and y directions 1529b1a0a8a3SJed Brown . dof - degrees of freedom per node 1530b1a0a8a3SJed Brown - overlap - overlap in mesh lines 1531b1a0a8a3SJed Brown 1532b1a0a8a3SJed Brown Output Parameters: 1533af538404SDmitry Karpeev + Nsub - the number of local subdomains created 15346a4f0f73SDmitry Karpeev . iis - array of index sets defining inner (nonoverlapping) subdomains 15356a4f0f73SDmitry Karpeev - ois - array of index sets defining outer (overlapping, if overlap > 0) subdomains 1536b1a0a8a3SJed Brown 1537b1a0a8a3SJed Brown 1538b1a0a8a3SJed Brown Level: advanced 1539b1a0a8a3SJed Brown 1540f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, create, subdomains, 2D, regular grid 1541b1a0a8a3SJed Brown 15428f3b4b4dSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMGetSubKSP(), PCGASMSetOverlap() 1543b1a0a8a3SJed Brown @*/ 15446a4f0f73SDmitry Karpeev PetscErrorCode PCGASMCreateSubdomains2D(PC pc,PetscInt M,PetscInt N,PetscInt Mdomains,PetscInt Ndomains,PetscInt dof,PetscInt overlap,PetscInt *nsub,IS **iis,IS **ois) 1545b1a0a8a3SJed Brown { 1546b1a0a8a3SJed Brown PetscErrorCode ierr; 15472bbb417fSDmitry Karpeev PetscMPIInt size, rank; 15482bbb417fSDmitry Karpeev PetscInt i, j; 15492bbb417fSDmitry Karpeev PetscInt maxheight, maxwidth; 15502bbb417fSDmitry Karpeev PetscInt xstart, xleft, xright, xleft_loc, xright_loc; 15512bbb417fSDmitry Karpeev PetscInt ystart, ylow, yhigh, ylow_loc, yhigh_loc; 1552eec7e3faSDmitry Karpeev PetscInt x[2][2], y[2][2], n[2]; 15532bbb417fSDmitry Karpeev PetscInt first, last; 15542bbb417fSDmitry Karpeev PetscInt nidx, *idx; 15552bbb417fSDmitry Karpeev PetscInt ii,jj,s,q,d; 1556af538404SDmitry Karpeev PetscInt k,kk; 15572bbb417fSDmitry Karpeev PetscMPIInt color; 15582bbb417fSDmitry Karpeev MPI_Comm comm, subcomm; 15596a4f0f73SDmitry Karpeev IS **xis = 0, **is = ois, **is_local = iis; 1560d34fcf5fSBarry Smith 1561b1a0a8a3SJed Brown PetscFunctionBegin; 15622bbb417fSDmitry Karpeev ierr = PetscObjectGetComm((PetscObject)pc, &comm);CHKERRQ(ierr); 15632bbb417fSDmitry Karpeev ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 15642bbb417fSDmitry Karpeev ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 15652bbb417fSDmitry Karpeev ierr = MatGetOwnershipRange(pc->pmat, &first, &last);CHKERRQ(ierr); 1566d34fcf5fSBarry Smith if (first%dof || last%dof) SETERRQ3(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Matrix row partitioning unsuitable for domain decomposition: local row range (%D,%D) " 15672bbb417fSDmitry Karpeev "does not respect the number of degrees of freedom per grid point %D", first, last, dof); 1568d34fcf5fSBarry Smith 1569af538404SDmitry Karpeev /* Determine the number of domains with nonzero intersections with the local ownership range. */ 15702bbb417fSDmitry Karpeev s = 0; 1571af538404SDmitry Karpeev ystart = 0; 1572af538404SDmitry Karpeev for (j=0; j<Ndomains; ++j) { 1573af538404SDmitry Karpeev maxheight = N/Ndomains + ((N % Ndomains) > j); /* Maximal height of subdomain */ 1574af538404SDmitry Karpeev if (maxheight < 2) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many %D subdomains in the vertical directon for mesh height %D", Ndomains, N); 1575eec7e3faSDmitry Karpeev /* Vertical domain limits with an overlap. */ 1576eec7e3faSDmitry Karpeev ylow = PetscMax(ystart - overlap,0); 1577eec7e3faSDmitry Karpeev yhigh = PetscMin(ystart + maxheight + overlap,N); 1578b1a0a8a3SJed Brown xstart = 0; 1579af538404SDmitry Karpeev for (i=0; i<Mdomains; ++i) { 1580af538404SDmitry Karpeev maxwidth = M/Mdomains + ((M % Mdomains) > i); /* Maximal width of subdomain */ 1581af538404SDmitry Karpeev if (maxwidth < 2) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many %D subdomains in the horizontal direction for mesh width %D", Mdomains, M); 1582eec7e3faSDmitry Karpeev /* Horizontal domain limits with an overlap. */ 1583eec7e3faSDmitry Karpeev xleft = PetscMax(xstart - overlap,0); 1584eec7e3faSDmitry Karpeev xright = PetscMin(xstart + maxwidth + overlap,M); 1585af538404SDmitry Karpeev /* 1586af538404SDmitry Karpeev Determine whether this subdomain intersects this processor's ownership range of pc->pmat. 1587af538404SDmitry Karpeev */ 1588c3518bceSDmitry Karpeev PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,(&xleft_loc),(&ylow_loc),(&xright_loc),(&yhigh_loc),(&nidx)); 15892fa5cd67SKarl Rupp if (nidx) ++s; 1590af538404SDmitry Karpeev xstart += maxwidth; 1591af538404SDmitry Karpeev } /* for (i = 0; i < Mdomains; ++i) */ 1592af538404SDmitry Karpeev ystart += maxheight; 1593af538404SDmitry Karpeev } /* for (j = 0; j < Ndomains; ++j) */ 15942fa5cd67SKarl Rupp 1595af538404SDmitry Karpeev /* Now we can allocate the necessary number of ISs. */ 1596af538404SDmitry Karpeev *nsub = s; 1597854ce69bSBarry Smith ierr = PetscMalloc1(*nsub,is);CHKERRQ(ierr); 1598854ce69bSBarry Smith ierr = PetscMalloc1(*nsub,is_local);CHKERRQ(ierr); 1599af538404SDmitry Karpeev s = 0; 1600af538404SDmitry Karpeev ystart = 0; 1601af538404SDmitry Karpeev for (j=0; j<Ndomains; ++j) { 1602af538404SDmitry Karpeev maxheight = N/Ndomains + ((N % Ndomains) > j); /* Maximal height of subdomain */ 1603af538404SDmitry Karpeev if (maxheight < 2) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many %D subdomains in the vertical directon for mesh height %D", Ndomains, N); 1604af538404SDmitry Karpeev /* Vertical domain limits with an overlap. */ 1605eec7e3faSDmitry Karpeev y[0][0] = PetscMax(ystart - overlap,0); 1606eec7e3faSDmitry Karpeev y[0][1] = PetscMin(ystart + maxheight + overlap,N); 1607eec7e3faSDmitry Karpeev /* Vertical domain limits without an overlap. */ 1608eec7e3faSDmitry Karpeev y[1][0] = ystart; 1609eec7e3faSDmitry Karpeev y[1][1] = PetscMin(ystart + maxheight,N); 1610eec7e3faSDmitry Karpeev xstart = 0; 1611af538404SDmitry Karpeev for (i=0; i<Mdomains; ++i) { 1612af538404SDmitry Karpeev maxwidth = M/Mdomains + ((M % Mdomains) > i); /* Maximal width of subdomain */ 1613af538404SDmitry Karpeev if (maxwidth < 2) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many %D subdomains in the horizontal direction for mesh width %D", Mdomains, M); 1614af538404SDmitry Karpeev /* Horizontal domain limits with an overlap. */ 1615eec7e3faSDmitry Karpeev x[0][0] = PetscMax(xstart - overlap,0); 1616eec7e3faSDmitry Karpeev x[0][1] = PetscMin(xstart + maxwidth + overlap,M); 1617eec7e3faSDmitry Karpeev /* Horizontal domain limits without an overlap. */ 1618eec7e3faSDmitry Karpeev x[1][0] = xstart; 1619eec7e3faSDmitry Karpeev x[1][1] = PetscMin(xstart+maxwidth,M); 16202bbb417fSDmitry Karpeev /* 16212bbb417fSDmitry Karpeev Determine whether this domain intersects this processor's ownership range of pc->pmat. 16222bbb417fSDmitry Karpeev Do this twice: first for the domains with overlaps, and once without. 16232bbb417fSDmitry Karpeev During the first pass create the subcommunicators, and use them on the second pass as well. 16242bbb417fSDmitry Karpeev */ 16252bbb417fSDmitry Karpeev for (q = 0; q < 2; ++q) { 1626cc96b2e5SBarry Smith PetscBool split = PETSC_FALSE; 16272bbb417fSDmitry Karpeev /* 16282bbb417fSDmitry Karpeev domain limits, (xleft, xright) and (ylow, yheigh) are adjusted 16292bbb417fSDmitry Karpeev according to whether the domain with an overlap or without is considered. 16302bbb417fSDmitry Karpeev */ 16312bbb417fSDmitry Karpeev xleft = x[q][0]; xright = x[q][1]; 16322bbb417fSDmitry Karpeev ylow = y[q][0]; yhigh = y[q][1]; 1633c3518bceSDmitry Karpeev PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,(&xleft_loc),(&ylow_loc),(&xright_loc),(&yhigh_loc),(&nidx)); 1634af538404SDmitry Karpeev nidx *= dof; 1635eec7e3faSDmitry Karpeev n[q] = nidx; 16362bbb417fSDmitry Karpeev /* 1637eec7e3faSDmitry Karpeev Based on the counted number of indices in the local domain *with an overlap*, 1638af538404SDmitry Karpeev construct a subcommunicator of all the processors supporting this domain. 1639eec7e3faSDmitry Karpeev Observe that a domain with an overlap might have nontrivial local support, 1640eec7e3faSDmitry Karpeev while the domain without an overlap might not. Hence, the decision to participate 1641eec7e3faSDmitry Karpeev in the subcommunicator must be based on the domain with an overlap. 16422bbb417fSDmitry Karpeev */ 16432bbb417fSDmitry Karpeev if (q == 0) { 16442fa5cd67SKarl Rupp if (nidx) color = 1; 16452fa5cd67SKarl Rupp else color = MPI_UNDEFINED; 16462bbb417fSDmitry Karpeev ierr = MPI_Comm_split(comm, color, rank, &subcomm);CHKERRQ(ierr); 1647cc96b2e5SBarry Smith split = PETSC_TRUE; 16482bbb417fSDmitry Karpeev } 1649af538404SDmitry Karpeev /* 1650eec7e3faSDmitry Karpeev Proceed only if the number of local indices *with an overlap* is nonzero. 1651af538404SDmitry Karpeev */ 1652eec7e3faSDmitry Karpeev if (n[0]) { 16532fa5cd67SKarl Rupp if (q == 0) xis = is; 1654af538404SDmitry Karpeev if (q == 1) { 1655af538404SDmitry Karpeev /* 1656eec7e3faSDmitry Karpeev The IS for the no-overlap subdomain shares a communicator with the overlapping domain. 1657af538404SDmitry Karpeev Moreover, if the overlap is zero, the two ISs are identical. 1658af538404SDmitry Karpeev */ 1659b1a0a8a3SJed Brown if (overlap == 0) { 1660eec7e3faSDmitry Karpeev (*is_local)[s] = (*is)[s]; 16612bbb417fSDmitry Karpeev ierr = PetscObjectReference((PetscObject)(*is)[s]);CHKERRQ(ierr); 16622bbb417fSDmitry Karpeev continue; 1663d34fcf5fSBarry Smith } else { 16646a4f0f73SDmitry Karpeev xis = is_local; 1665eec7e3faSDmitry Karpeev subcomm = ((PetscObject)(*is)[s])->comm; 16662bbb417fSDmitry Karpeev } 1667af538404SDmitry Karpeev } /* if (q == 1) */ 16680298fd71SBarry Smith idx = NULL; 1669785e854fSJed Brown ierr = PetscMalloc1(nidx,&idx);CHKERRQ(ierr); 1670eec7e3faSDmitry Karpeev if (nidx) { 16712bbb417fSDmitry Karpeev k = 0; 16722bbb417fSDmitry Karpeev for (jj=ylow_loc; jj<yhigh_loc; ++jj) { 1673af538404SDmitry Karpeev PetscInt x0 = (jj==ylow_loc) ? xleft_loc : xleft; 1674af538404SDmitry Karpeev PetscInt x1 = (jj==yhigh_loc-1) ? xright_loc : xright; 1675af538404SDmitry Karpeev kk = dof*(M*jj + x0); 16762bbb417fSDmitry Karpeev for (ii=x0; ii<x1; ++ii) { 16772bbb417fSDmitry Karpeev for (d = 0; d < dof; ++d) { 16782bbb417fSDmitry Karpeev idx[k++] = kk++; 1679b1a0a8a3SJed Brown } 1680b1a0a8a3SJed Brown } 1681b1a0a8a3SJed Brown } 1682eec7e3faSDmitry Karpeev } 16836a4f0f73SDmitry Karpeev ierr = ISCreateGeneral(subcomm,nidx,idx,PETSC_OWN_POINTER,(*xis)+s);CHKERRQ(ierr); 1684cc96b2e5SBarry Smith if (split) { 1685cc96b2e5SBarry Smith ierr = MPI_Comm_free(&subcomm);CHKERRQ(ierr); 1686cc96b2e5SBarry Smith } 1687eec7e3faSDmitry Karpeev }/* if (n[0]) */ 16882bbb417fSDmitry Karpeev }/* for (q = 0; q < 2; ++q) */ 16892fa5cd67SKarl Rupp if (n[0]) ++s; 1690af538404SDmitry Karpeev xstart += maxwidth; 1691af538404SDmitry Karpeev } /* for (i = 0; i < Mdomains; ++i) */ 16922bbb417fSDmitry Karpeev ystart += maxheight; 1693af538404SDmitry Karpeev } /* for (j = 0; j < Ndomains; ++j) */ 1694b1a0a8a3SJed Brown PetscFunctionReturn(0); 1695b1a0a8a3SJed Brown } 1696b1a0a8a3SJed Brown 1697b1a0a8a3SJed Brown #undef __FUNCT__ 169806b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMGetSubdomains" 1699b1a0a8a3SJed Brown /*@C 170006b43e7eSDmitry Karpeev PCGASMGetSubdomains - Gets the subdomains supported on this processor 170106b43e7eSDmitry Karpeev for the additive Schwarz preconditioner. 1702b1a0a8a3SJed Brown 1703b1a0a8a3SJed Brown Not Collective 1704b1a0a8a3SJed Brown 1705b1a0a8a3SJed Brown Input Parameter: 1706b1a0a8a3SJed Brown . pc - the preconditioner context 1707b1a0a8a3SJed Brown 1708b1a0a8a3SJed Brown Output Parameters: 1709b1a0a8a3SJed Brown + n - the number of subdomains for this processor (default value = 1) 17100298fd71SBarry Smith . iis - the index sets that define the inner subdomains (without overlap) supported on this processor (can be NULL) 17110298fd71SBarry Smith - ois - the index sets that define the outer subdomains (with overlap) supported on this processor (can be NULL) 1712b1a0a8a3SJed Brown 1713b1a0a8a3SJed Brown 1714b1a0a8a3SJed Brown Notes: 17156a4f0f73SDmitry Karpeev The user is responsible for destroying the ISs and freeing the returned arrays. 1716b1a0a8a3SJed Brown The IS numbering is in the parallel, global numbering of the vector. 1717b1a0a8a3SJed Brown 1718b1a0a8a3SJed Brown Level: advanced 1719b1a0a8a3SJed Brown 172006b43e7eSDmitry Karpeev .keywords: PC, GASM, get, subdomains, additive Schwarz 1721b1a0a8a3SJed Brown 17228f3b4b4dSDmitry Karpeev .seealso: PCGASMSetOverlap(), PCGASMGetSubKSP(), PCGASMCreateSubdomains2D(), 17238f3b4b4dSDmitry Karpeev PCGASMSetSubdomains(), PCGASMGetSubmatrices() 1724b1a0a8a3SJed Brown @*/ 17256a4f0f73SDmitry Karpeev PetscErrorCode PCGASMGetSubdomains(PC pc,PetscInt *n,IS *iis[],IS *ois[]) 1726b1a0a8a3SJed Brown { 1727f746d493SDmitry Karpeev PC_GASM *osm; 1728b1a0a8a3SJed Brown PetscErrorCode ierr; 1729b1a0a8a3SJed Brown PetscBool match; 17306a4f0f73SDmitry Karpeev PetscInt i; 17315fd66863SKarl Rupp 1732b1a0a8a3SJed Brown PetscFunctionBegin; 1733b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1734251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match);CHKERRQ(ierr); 1735ce94432eSBarry Smith if (!match) SETERRQ2(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Incorrect object type: expected %s, got %s instead", PCGASM, ((PetscObject)pc)->type_name); 1736f746d493SDmitry Karpeev osm = (PC_GASM*)pc->data; 1737ab3e923fSDmitry Karpeev if (n) *n = osm->n; 17386a4f0f73SDmitry Karpeev if (iis) { 1739785e854fSJed Brown ierr = PetscMalloc1(osm->n, iis);CHKERRQ(ierr); 17406a4f0f73SDmitry Karpeev } 17416a4f0f73SDmitry Karpeev if (ois) { 1742785e854fSJed Brown ierr = PetscMalloc1(osm->n, ois);CHKERRQ(ierr); 17436a4f0f73SDmitry Karpeev } 17446a4f0f73SDmitry Karpeev if (iis || ois) { 17456a4f0f73SDmitry Karpeev for (i = 0; i < osm->n; ++i) { 17466a4f0f73SDmitry Karpeev if (iis) (*iis)[i] = osm->iis[i]; 17476a4f0f73SDmitry Karpeev if (ois) (*ois)[i] = osm->ois[i]; 17486a4f0f73SDmitry Karpeev } 1749b1a0a8a3SJed Brown } 1750b1a0a8a3SJed Brown PetscFunctionReturn(0); 1751b1a0a8a3SJed Brown } 1752b1a0a8a3SJed Brown 1753b1a0a8a3SJed Brown #undef __FUNCT__ 175406b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMGetSubmatrices" 1755b1a0a8a3SJed Brown /*@C 175606b43e7eSDmitry Karpeev PCGASMGetSubmatrices - Gets the local submatrices (for this processor 1757b1a0a8a3SJed Brown only) for the additive Schwarz preconditioner. 1758b1a0a8a3SJed Brown 1759b1a0a8a3SJed Brown Not Collective 1760b1a0a8a3SJed Brown 1761b1a0a8a3SJed Brown Input Parameter: 1762b1a0a8a3SJed Brown . pc - the preconditioner context 1763b1a0a8a3SJed Brown 1764b1a0a8a3SJed Brown Output Parameters: 1765b1a0a8a3SJed Brown + n - the number of matrices for this processor (default value = 1) 1766b1a0a8a3SJed Brown - mat - the matrices 1767b1a0a8a3SJed Brown 176806b43e7eSDmitry Karpeev Notes: matrices returned by this routine have the same communicators as the index sets (IS) 17698f3b4b4dSDmitry Karpeev used to define subdomains in PCGASMSetSubdomains() 1770b1a0a8a3SJed Brown Level: advanced 1771b1a0a8a3SJed Brown 1772f746d493SDmitry Karpeev .keywords: PC, GASM, set, local, subdomains, additive Schwarz, block Jacobi 1773b1a0a8a3SJed Brown 17748f3b4b4dSDmitry Karpeev .seealso: PCGASMSetOverlap(), PCGASMGetSubKSP(), 177506b43e7eSDmitry Karpeev PCGASMCreateSubdomains2D(), PCGASMSetSubdomains(), PCGASMGetSubdomains() 1776b1a0a8a3SJed Brown @*/ 177706b43e7eSDmitry Karpeev PetscErrorCode PCGASMGetSubmatrices(PC pc,PetscInt *n,Mat *mat[]) 1778b1a0a8a3SJed Brown { 1779f746d493SDmitry Karpeev PC_GASM *osm; 1780b1a0a8a3SJed Brown PetscErrorCode ierr; 1781b1a0a8a3SJed Brown PetscBool match; 1782b1a0a8a3SJed Brown 1783b1a0a8a3SJed Brown PetscFunctionBegin; 1784b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1785b1a0a8a3SJed Brown PetscValidIntPointer(n,2); 1786b1a0a8a3SJed Brown if (mat) PetscValidPointer(mat,3); 1787ce94432eSBarry Smith if (!pc->setupcalled) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call after KSPSetUP() or PCSetUp()."); 1788251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match);CHKERRQ(ierr); 1789ce94432eSBarry Smith if (!match) SETERRQ2(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Expected %s, got %s instead", PCGASM, ((PetscObject)pc)->type_name); 1790f746d493SDmitry Karpeev osm = (PC_GASM*)pc->data; 1791ab3e923fSDmitry Karpeev if (n) *n = osm->n; 1792b1a0a8a3SJed Brown if (mat) *mat = osm->pmat; 1793b1a0a8a3SJed Brown PetscFunctionReturn(0); 1794b1a0a8a3SJed Brown } 1795d709ea83SDmitry Karpeev 1796d709ea83SDmitry Karpeev #undef __FUNCT__ 17978f3b4b4dSDmitry Karpeev #define __FUNCT__ "PCGASMSetUseDMSubdomains" 1798d709ea83SDmitry Karpeev /*@ 17998f3b4b4dSDmitry Karpeev PCGASMSetUseDMSubdomains - Indicates whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible. 1800d709ea83SDmitry Karpeev Logically Collective 1801d709ea83SDmitry Karpeev 1802d709ea83SDmitry Karpeev Input Parameter: 1803d709ea83SDmitry Karpeev + pc - the preconditioner 1804d709ea83SDmitry Karpeev - flg - boolean indicating whether to use subdomains defined by the DM 1805d709ea83SDmitry Karpeev 1806d709ea83SDmitry Karpeev Options Database Key: 18078f3b4b4dSDmitry Karpeev . -pc_gasm_dm_subdomains -pc_gasm_overlap -pc_gasm_total_subdomains 1808d709ea83SDmitry Karpeev 1809d709ea83SDmitry Karpeev Level: intermediate 1810d709ea83SDmitry Karpeev 1811d709ea83SDmitry Karpeev Notes: 18128f3b4b4dSDmitry Karpeev PCGASMSetSubdomains(), PCGASMSetTotalSubdomains() or PCGASMSetOverlap() take precedence over PCGASMSetUseDMSubdomains(), 18138f3b4b4dSDmitry Karpeev so setting PCGASMSetSubdomains() with nontrivial subdomain ISs or any of PCGASMSetTotalSubdomains() and PCGASMSetOverlap() 18148f3b4b4dSDmitry Karpeev automatically turns the latter off. 1815d709ea83SDmitry Karpeev 1816d709ea83SDmitry Karpeev .keywords: PC, ASM, DM, set, subdomains, additive Schwarz 1817d709ea83SDmitry Karpeev 18188f3b4b4dSDmitry Karpeev .seealso: PCGASMGetUseDMSubdomains(), PCGASMSetSubdomains(), PCGASMSetOverlap() 1819d709ea83SDmitry Karpeev PCGASMCreateSubdomains2D() 1820d709ea83SDmitry Karpeev @*/ 18218f3b4b4dSDmitry Karpeev PetscErrorCode PCGASMSetUseDMSubdomains(PC pc,PetscBool flg) 1822d709ea83SDmitry Karpeev { 1823d709ea83SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 1824d709ea83SDmitry Karpeev PetscErrorCode ierr; 1825d709ea83SDmitry Karpeev PetscBool match; 1826d709ea83SDmitry Karpeev 1827d709ea83SDmitry Karpeev PetscFunctionBegin; 1828d709ea83SDmitry Karpeev PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1829d709ea83SDmitry Karpeev PetscValidLogicalCollectiveBool(pc,flg,2); 1830d709ea83SDmitry Karpeev if (pc->setupcalled) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Not for a setup PC."); 1831d709ea83SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match);CHKERRQ(ierr); 1832d709ea83SDmitry Karpeev if (match) { 18338f3b4b4dSDmitry Karpeev if (!osm->user_subdomains && osm->N == PETSC_DETERMINE && osm->overlap < 0) { 1834d709ea83SDmitry Karpeev osm->dm_subdomains = flg; 1835d709ea83SDmitry Karpeev } 18368f3b4b4dSDmitry Karpeev } 1837d709ea83SDmitry Karpeev PetscFunctionReturn(0); 1838d709ea83SDmitry Karpeev } 1839d709ea83SDmitry Karpeev 1840d709ea83SDmitry Karpeev #undef __FUNCT__ 18418f3b4b4dSDmitry Karpeev #define __FUNCT__ "PCGASMGetUseDMSubdomains" 1842d709ea83SDmitry Karpeev /*@ 18438f3b4b4dSDmitry Karpeev PCGASMGetUseDMSubdomains - Returns flag indicating whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible. 1844d709ea83SDmitry Karpeev Not Collective 1845d709ea83SDmitry Karpeev 1846d709ea83SDmitry Karpeev Input Parameter: 1847d709ea83SDmitry Karpeev . pc - the preconditioner 1848d709ea83SDmitry Karpeev 1849d709ea83SDmitry Karpeev Output Parameter: 1850d709ea83SDmitry Karpeev . flg - boolean indicating whether to use subdomains defined by the DM 1851d709ea83SDmitry Karpeev 1852d709ea83SDmitry Karpeev Level: intermediate 1853d709ea83SDmitry Karpeev 1854d709ea83SDmitry Karpeev .keywords: PC, ASM, DM, set, subdomains, additive Schwarz 1855d709ea83SDmitry Karpeev 18568f3b4b4dSDmitry Karpeev .seealso: PCGASMSetUseDMSubdomains(), PCGASMSetOverlap() 1857d709ea83SDmitry Karpeev PCGASMCreateSubdomains2D() 1858d709ea83SDmitry Karpeev @*/ 18598f3b4b4dSDmitry Karpeev PetscErrorCode PCGASMGetUseDMSubdomains(PC pc,PetscBool* flg) 1860d709ea83SDmitry Karpeev { 1861d709ea83SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 1862d709ea83SDmitry Karpeev PetscErrorCode ierr; 1863d709ea83SDmitry Karpeev PetscBool match; 1864d709ea83SDmitry Karpeev 1865d709ea83SDmitry Karpeev PetscFunctionBegin; 1866d709ea83SDmitry Karpeev PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1867d709ea83SDmitry Karpeev PetscValidPointer(flg,2); 1868d709ea83SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match);CHKERRQ(ierr); 1869d709ea83SDmitry Karpeev if (match) { 1870d709ea83SDmitry Karpeev if (flg) *flg = osm->dm_subdomains; 1871d709ea83SDmitry Karpeev } 1872d709ea83SDmitry Karpeev PetscFunctionReturn(0); 1873d709ea83SDmitry Karpeev } 1874