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 7*8f3b4b4dSDmitry 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()) 9*8f3b4b4dSDmitry 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 */ 17*8f3b4b4dSDmitry Karpeev PCGASMType type; /* use reduced interpolation, restriction or both */ 18*8f3b4b4dSDmitry Karpeev PetscBool type_set; /* if user set this value (so won't change it for symmetric problems) */ 19*8f3b4b4dSDmitry Karpeev PetscBool same_subdomain_solvers; /* flag indicating whether all local solvers are same */ 20*8f3b4b4dSDmitry Karpeev PetscBool sort_indices; /* flag to sort subdomain indices */ 21*8f3b4b4dSDmitry Karpeev PetscBool user_subdomains; /* whether the user set explicit subdomain index sets -- keep them on PCReset() */ 22*8f3b4b4dSDmitry Karpeev PetscBool dm_subdomains; /* whether DM is allowed to define subdomains */ 23*8f3b4b4dSDmitry Karpeev IS *ois; /* index sets that define the outer (conceptually, overlapping) subdomains */ 24*8f3b4b4dSDmitry Karpeev IS *iis; /* index sets that define the inner (conceptually, nonoverlapping) subdomains */ 25*8f3b4b4dSDmitry Karpeev KSP *ksp; /* linear solvers for each subdomain */ 26*8f3b4b4dSDmitry 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__ 34*8f3b4b4dSDmitry Karpeev #define __FUNCT__ "PCGASMComputeGlobalSubdomainNumbering_Private" 35*8f3b4b4dSDmitry Karpeev static PetscErrorCode PCGASMComputeGlobalSubdomainNumbering_Private(PC pc,PetscInt **numbering,PetscInt **permutation) 36*8f3b4b4dSDmitry Karpeev { 37*8f3b4b4dSDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 38*8f3b4b4dSDmitry Karpeev PetscInt i; 39*8f3b4b4dSDmitry Karpeev PetscErrorCode ierr; 40*8f3b4b4dSDmitry Karpeev 41*8f3b4b4dSDmitry Karpeev PetscFunctionBegin; 42*8f3b4b4dSDmitry Karpeev /* Determine the number of globally-distinct subdomains and compute a global numbering for them. */ 43*8f3b4b4dSDmitry Karpeev ierr = PetscMalloc2(osm->n,numbering,osm->n,permutation);CHKERRQ(ierr); 44*8f3b4b4dSDmitry Karpeev ierr = PetscObjectsListGetGlobalNumbering(PetscObjectComm((PetscObject)pc),osm->n,(PetscObject*)osm->iis,NULL,*numbering);CHKERRQ(ierr); 45*8f3b4b4dSDmitry Karpeev for (i = 0; i < osm->n; ++i) (*permutation)[i] = i; 46*8f3b4b4dSDmitry Karpeev ierr = PetscSortIntWithPermutation(osm->n,*numbering,*permutation);CHKERRQ(ierr); 47*8f3b4b4dSDmitry Karpeev PetscFunctionReturn(0); 48*8f3b4b4dSDmitry Karpeev } 49*8f3b4b4dSDmitry Karpeev 50*8f3b4b4dSDmitry 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); 814bde246dSDmitry Karpeev ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr); 824bde246dSDmitry Karpeev ierr = PetscViewerASCIISynchronizedPrintf(viewer, "%s", cidx);CHKERRQ(ierr); 834bde246dSDmitry Karpeev ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 844bde246dSDmitry Karpeev ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);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); 1067b23a99aSBarry Smith ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr); 107af538404SDmitry Karpeev ierr = PetscViewerASCIISynchronizedPrintf(viewer, "%s", cidx);CHKERRQ(ierr); 108af538404SDmitry Karpeev ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 1097b23a99aSBarry Smith ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);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]; 123*8f3b4b4dSDmitry Karpeev PetscInt l, d, count; 124*8f3b4b4dSDmitry Karpeev PetscBool doprint,found; 1250298fd71SBarry Smith PetscViewer viewer, sviewer = NULL; 126*8f3b4b4dSDmitry 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); 131*8f3b4b4dSDmitry Karpeev doprint = PETSC_FALSE; 132*8f3b4b4dSDmitry Karpeev ierr = PetscOptionsGetBool(prefix,"-pc_gasm_print_subdomains",&doprint,NULL);CHKERRQ(ierr); 133*8f3b4b4dSDmitry 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; 143*8f3b4b4dSDmitry Karpeev ierr = PCGASMComputeGlobalSubdomainNumbering_Private(pc,&numbering,&permutation);CHKERRQ(ierr); 144*8f3b4b4dSDmitry 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) { 1494bde246dSDmitry Karpeev ierr = PetscViewerGetSubcomm(viewer,((PetscObject)osm->ois[d])->comm, &sviewer);CHKERRQ(ierr); 1504bde246dSDmitry Karpeev ierr = PCGASMSubdomainView_Private(pc,d,sviewer);CHKERRQ(ierr); 1514bde246dSDmitry Karpeev ierr = PetscViewerRestoreSubcomm(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 } 157*8f3b4b4dSDmitry 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; 171*8f3b4b4dSDmitry Karpeev PetscInt bsz; 17206b43e7eSDmitry Karpeev PetscBool iascii,view_subdomains=PETSC_FALSE; 173b1a0a8a3SJed Brown PetscViewer sviewer; 174*8f3b4b4dSDmitry 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"; 178*8f3b4b4dSDmitry 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 } 187*8f3b4b4dSDmitry Karpeev if (osm->N != PETSC_DETERMINE) { 188*8f3b4b4dSDmitry Karpeev ierr = PetscSNPrintf(gsubdomains, sizeof(gsubdomains), "total number of subdomains = %D",osm->N);CHKERRQ(ierr); 18906b43e7eSDmitry Karpeev } 190*8f3b4b4dSDmitry 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); 2107b23a99aSBarry Smith ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr); 211*8f3b4b4dSDmitry 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); 2137b23a99aSBarry Smith ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);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. */ 221*8f3b4b4dSDmitry Karpeev ierr = PCGASMComputeGlobalSubdomainNumbering_Private(pc,&numbering,&permutation);CHKERRQ(ierr); 22206b43e7eSDmitry Karpeev l = 0; 223*8f3b4b4dSDmitry 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); 2306a4f0f73SDmitry Karpeev ierr = PetscViewerGetSubcomm(viewer,((PetscObject)osm->ois[d])->comm, &sviewer);CHKERRQ(ierr); 2316a4f0f73SDmitry Karpeev ierr = ISGetLocalSize(osm->ois[d],&bsz);CHKERRQ(ierr); 2327b23a99aSBarry Smith ierr = PetscViewerASCIISynchronizedAllow(sviewer,PETSC_TRUE);CHKERRQ(ierr); 233*8f3b4b4dSDmitry 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); 2346a4f0f73SDmitry Karpeev ierr = PetscViewerFlush(sviewer);CHKERRQ(ierr); 2356a4f0f73SDmitry Karpeev ierr = PetscViewerASCIISynchronizedAllow(sviewer,PETSC_FALSE);CHKERRQ(ierr); 23606b43e7eSDmitry Karpeev if (view_subdomains) { 23706b43e7eSDmitry Karpeev ierr = PCGASMSubdomainView_Private(pc,d,sviewer);CHKERRQ(ierr); 23806b43e7eSDmitry Karpeev } 2396a4f0f73SDmitry Karpeev if (!pc->setupcalled) { 2406a4f0f73SDmitry Karpeev PetscViewerASCIIPrintf(sviewer, "Solver not set up yet: PCSetUp() not yet called\n");CHKERRQ(ierr); 2412fa5cd67SKarl Rupp } else { 24206b43e7eSDmitry Karpeev ierr = KSPView(osm->ksp[d],sviewer);CHKERRQ(ierr); 2436a4f0f73SDmitry Karpeev } 24406b43e7eSDmitry Karpeev ierr = PetscViewerASCIIPrintf(sviewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr); 2457b23a99aSBarry Smith ierr = PetscViewerFlush(sviewer);CHKERRQ(ierr); 2466a4f0f73SDmitry Karpeev ierr = PetscViewerRestoreSubcomm(viewer,((PetscObject)osm->ois[d])->comm, &sviewer);CHKERRQ(ierr); 24706b43e7eSDmitry Karpeev ++l; 248b1a0a8a3SJed Brown } 24906b43e7eSDmitry Karpeev } 250ce94432eSBarry Smith ierr = MPI_Barrier(PetscObjectComm((PetscObject)pc));CHKERRQ(ierr); 251b1a0a8a3SJed Brown } 252*8f3b4b4dSDmitry Karpeev ierr = PetscFree2(numbering,permutation);CHKERRQ(ierr); 253b1a0a8a3SJed Brown ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 254b1a0a8a3SJed Brown } 255b1a0a8a3SJed Brown PetscFunctionReturn(0); 256b1a0a8a3SJed Brown } 257b1a0a8a3SJed Brown 258*8f3b4b4dSDmitry 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; 2812fa5cd67SKarl Rupp PetscInt gofirst; /* Start of locally-owned indices in the vectors -- osm->gx,osm->gy -- over the disjoint union of outer subdomains. */ 282*8f3b4b4dSDmitry Karpeev PetscInt num_subdomains = 0; 2830298fd71SBarry Smith DM *subdomain_dm = NULL; 284*8f3b4b4dSDmitry Karpeev char **subdomain_names = NULL; 285*8f3b4b4dSDmitry Karpeev 286b1a0a8a3SJed Brown 287b1a0a8a3SJed Brown PetscFunctionBegin; 288ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr); 289ce94432eSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);CHKERRQ(ierr); 290b1a0a8a3SJed Brown if (!pc->setupcalled) { 291b1a0a8a3SJed Brown 292b1a0a8a3SJed Brown if (!osm->type_set) { 293b1a0a8a3SJed Brown ierr = MatIsSymmetricKnown(pc->pmat,&symset,&flg);CHKERRQ(ierr); 2942fa5cd67SKarl Rupp if (symset && flg) osm->type = PC_GASM_BASIC; 295b1a0a8a3SJed Brown } 296b1a0a8a3SJed Brown 297*8f3b4b4dSDmitry Karpeev if (osm->n == PETSC_DETERMINE) { 298*8f3b4b4dSDmitry Karpeev if (osm->N != PETSC_DETERMINE) { 299*8f3b4b4dSDmitry Karpeev /* No local subdomains given, but the desired number of total subdomains is known, so construct them accordingly. */ 300*8f3b4b4dSDmitry Karpeev ierr = PCGASMCreateSubdomains(pc->pmat,osm->N,&osm->n,&osm->iis);CHKERRQ(ierr); 301*8f3b4b4dSDmitry Karpeev } else if (osm->dm_subdomains && pc->dm) { 302*8f3b4b4dSDmitry Karpeev /* try pc->dm next, if allowed */ 303*8f3b4b4dSDmitry Karpeev PetscInt d; 304a35b7b57SDmitry Karpeev IS *inner_subdomain_is, *outer_subdomain_is; 305a35b7b57SDmitry Karpeev ierr = DMCreateDomainDecomposition(pc->dm, &num_subdomains, &subdomain_names, &inner_subdomain_is, &outer_subdomain_is, &subdomain_dm);CHKERRQ(ierr); 306a35b7b57SDmitry Karpeev if (num_subdomains) { 307a35b7b57SDmitry Karpeev ierr = PCGASMSetSubdomains(pc, num_subdomains, inner_subdomain_is, outer_subdomain_is);CHKERRQ(ierr); 30869ca1f37SDmitry Karpeev } 309a35b7b57SDmitry Karpeev for (d = 0; d < num_subdomains; ++d) { 310a35b7b57SDmitry Karpeev if (inner_subdomain_is) {ierr = ISDestroy(&inner_subdomain_is[d]);CHKERRQ(ierr);} 311a35b7b57SDmitry Karpeev if (outer_subdomain_is) {ierr = ISDestroy(&outer_subdomain_is[d]);CHKERRQ(ierr);} 312fdc48646SDmitry Karpeev } 313a35b7b57SDmitry Karpeev ierr = PetscFree(inner_subdomain_is);CHKERRQ(ierr); 314a35b7b57SDmitry Karpeev ierr = PetscFree(outer_subdomain_is);CHKERRQ(ierr); 315*8f3b4b4dSDmitry Karpeev } else { 316*8f3b4b4dSDmitry Karpeev /* still no subdomains; use one per processor */ 31744fe18e5SDmitry Karpeev osm->nmax = osm->n = 1; 318ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr); 319f746d493SDmitry Karpeev osm->N = size; 320*8f3b4b4dSDmitry Karpeev ierr = PCGASMCreateLocalSubdomains(pc->pmat,osm->n,&osm->iis);CHKERRQ(ierr); 321fdc48646SDmitry Karpeev } 32206b43e7eSDmitry Karpeev } 323a35b7b57SDmitry Karpeev if (!osm->iis) { 324a35b7b57SDmitry Karpeev /* 325*8f3b4b4dSDmitry Karpeev osm->n was set in PCGASMSetSubdomains(), but the actual subdomains have not been supplied. 326*8f3b4b4dSDmitry Karpeev We create the requisite number of local inner subdomains and then expand them into 327*8f3b4b4dSDmitry Karpeev out subdomains, if necessary. 328a35b7b57SDmitry Karpeev */ 329*8f3b4b4dSDmitry Karpeev ierr = PCGASMCreateLocalSubdomains(pc->pmat,osm->n,&osm->iis);CHKERRQ(ierr); 330a35b7b57SDmitry Karpeev } 331*8f3b4b4dSDmitry Karpeev if (!osm->ois) { 332*8f3b4b4dSDmitry Karpeev /* 333*8f3b4b4dSDmitry Karpeev Initially make outer subdomains the same as inner subdomains. If nonzero additional overlap 334*8f3b4b4dSDmitry Karpeev has been requested, copy the inner subdomains over so they can be modified. 335*8f3b4b4dSDmitry Karpeev */ 336*8f3b4b4dSDmitry Karpeev ierr = PetscMalloc1(osm->n,&osm->ois);CHKERRQ(ierr); 337*8f3b4b4dSDmitry Karpeev for (i=0; i<osm->n; ++i) { 338*8f3b4b4dSDmitry Karpeev if (osm->overlap > 0) { /* With positive overlap, osm->iis[i] will be modified */ 339*8f3b4b4dSDmitry Karpeev ierr = ISDuplicate(osm->iis[i],(osm->ois)+i);CHKERRQ(ierr); 340*8f3b4b4dSDmitry Karpeev ierr = ISCopy(osm->iis[i],osm->ois[i]);CHKERRQ(ierr); 341*8f3b4b4dSDmitry Karpeev } else { 342*8f3b4b4dSDmitry Karpeev ierr = PetscObjectReference((PetscObject)((osm->iis)[i]));CHKERRQ(ierr); 343*8f3b4b4dSDmitry Karpeev osm->ois[i] = osm->iis[i]; 344*8f3b4b4dSDmitry Karpeev } 345*8f3b4b4dSDmitry Karpeev } 346*8f3b4b4dSDmitry Karpeev if (osm->overlap > 0) { 347*8f3b4b4dSDmitry Karpeev /* Extend the "overlapping" regions by a number of steps */ 348*8f3b4b4dSDmitry Karpeev ierr = MatIncreaseOverlap(pc->pmat,osm->n,osm->ois,osm->overlap);CHKERRQ(ierr); 349*8f3b4b4dSDmitry Karpeev } 350b1a0a8a3SJed Brown } 3516a4f0f73SDmitry Karpeev 352*8f3b4b4dSDmitry Karpeev /* Now the subdomains are defined. Determine their global and max local numbers, if necessary. */ 353*8f3b4b4dSDmitry Karpeev if (osm->nmax == PETSC_DETERMINE) { 354*8f3b4b4dSDmitry Karpeev PetscMPIInt inwork,outwork; 355*8f3b4b4dSDmitry Karpeev /* determine global number of subdomains and the max number of local subdomains */ 356*8f3b4b4dSDmitry Karpeev inwork = osm->n; 357*8f3b4b4dSDmitry Karpeev ierr = MPI_Allreduce(&inwork,&outwork,1,MPI_INT,MPI_MAX,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr); 358*8f3b4b4dSDmitry Karpeev osm->nmax = outwork; 359*8f3b4b4dSDmitry Karpeev } 360*8f3b4b4dSDmitry Karpeev if (osm->N == PETSC_DETERMINE) { 361*8f3b4b4dSDmitry Karpeev /* Determine the number of globally-distinct subdomains and compute a global numbering for them. */ 362*8f3b4b4dSDmitry Karpeev ierr = PetscObjectsListGetGlobalNumbering(PetscObjectComm((PetscObject)pc),osm->n,(PetscObject*)osm->ois,&osm->N,NULL);CHKERRQ(ierr); 363*8f3b4b4dSDmitry Karpeev } 364*8f3b4b4dSDmitry Karpeev 365b1a0a8a3SJed Brown 366b1a0a8a3SJed Brown if (osm->sort_indices) { 367f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 3686a4f0f73SDmitry Karpeev ierr = ISSort(osm->ois[i]);CHKERRQ(ierr); 3696a4f0f73SDmitry Karpeev ierr = ISSort(osm->iis[i]);CHKERRQ(ierr); 370b1a0a8a3SJed Brown } 371b1a0a8a3SJed Brown } 372*8f3b4b4dSDmitry Karpeev 373*8f3b4b4dSDmitry Karpeev ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 374*8f3b4b4dSDmitry Karpeev ierr = PCGASMPrintSubdomains(pc);CHKERRQ(ierr); 375*8f3b4b4dSDmitry 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; 387f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 3886a4f0f73SDmitry Karpeev ierr = ISGetLocalSize(osm->ois[i],&oni);CHKERRQ(ierr); 3896a4f0f73SDmitry Karpeev ierr = ISGetIndices(osm->ois[i],&oidxi);CHKERRQ(ierr); 3906a4f0f73SDmitry Karpeev ierr = PetscMemcpy(oidx+on, oidxi, sizeof(PetscInt)*oni);CHKERRQ(ierr); 3916a4f0f73SDmitry Karpeev ierr = ISRestoreIndices(osm->ois[i], &oidxi);CHKERRQ(ierr); 3926a4f0f73SDmitry Karpeev on += oni; 393f746d493SDmitry Karpeev } 3946a4f0f73SDmitry Karpeev ierr = ISCreateGeneral(((PetscObject)(pc))->comm, on, oidx, PETSC_OWN_POINTER, &gois);CHKERRQ(ierr); 3952a7a6963SBarry Smith ierr = MatCreateVecs(pc->pmat,&x,&y);CHKERRQ(ierr); 396ce94432eSBarry Smith ierr = VecCreateMPI(PetscObjectComm((PetscObject)pc), on, PETSC_DECIDE, &osm->gx);CHKERRQ(ierr); 397f746d493SDmitry Karpeev ierr = VecDuplicate(osm->gx,&osm->gy);CHKERRQ(ierr); 3980298fd71SBarry Smith ierr = VecGetOwnershipRange(osm->gx, &gofirst, NULL);CHKERRQ(ierr); 399ce94432eSBarry Smith ierr = ISCreateStride(PetscObjectComm((PetscObject)pc),on,gofirst,1, &goid);CHKERRQ(ierr); 4006a4f0f73SDmitry Karpeev ierr = VecScatterCreate(x,gois,osm->gx,goid, &(osm->gorestriction));CHKERRQ(ierr); 4016a4f0f73SDmitry Karpeev ierr = VecDestroy(&x);CHKERRQ(ierr); 4027105d80bSDmitry Karpeev ierr = ISDestroy(&gois);CHKERRQ(ierr); 4032fa5cd67SKarl Rupp 4046a4f0f73SDmitry Karpeev /* Merge inner subdomain ISs and construct a restriction onto the disjoint union of local inner subdomains. */ 4052fa5cd67SKarl Rupp { 4062fa5cd67SKarl Rupp PetscInt ini; /* Number of indices the i-th a local inner subdomain. */ 4076a4f0f73SDmitry Karpeev PetscInt in; /* Number of indices in the disjoint uniont of local inner subdomains. */ 4086a4f0f73SDmitry Karpeev PetscInt *iidx; /* Global indices in the merged local inner subdomain. */ 4096a4f0f73SDmitry Karpeev PetscInt *ioidx; /* Global indices of the disjoint union of inner subdomains within the disjoint union of outer subdomains. */ 4106a4f0f73SDmitry Karpeev IS giis; /* IS for the disjoint union of inner subdomains. */ 4116a4f0f73SDmitry Karpeev IS giois; /* IS for the disjoint union of inner subdomains within the disjoint union of outer subdomains. */ 412f746d493SDmitry Karpeev /**/ 4136a4f0f73SDmitry Karpeev in = 0; 414f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 4156a4f0f73SDmitry Karpeev ierr = ISGetLocalSize(osm->iis[i],&ini);CHKERRQ(ierr); 4166a4f0f73SDmitry Karpeev in += ini; 417f746d493SDmitry Karpeev } 418785e854fSJed Brown ierr = PetscMalloc1(in, &iidx);CHKERRQ(ierr); 419785e854fSJed Brown ierr = PetscMalloc1(in, &ioidx);CHKERRQ(ierr); 4200298fd71SBarry Smith ierr = VecGetOwnershipRange(osm->gx,&gofirst, NULL);CHKERRQ(ierr); 4216a4f0f73SDmitry Karpeev in = 0; 4226a4f0f73SDmitry Karpeev on = 0; 423f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 4246a4f0f73SDmitry Karpeev const PetscInt *iidxi; /* Global indices of the i-th local inner subdomain. */ 4256a4f0f73SDmitry Karpeev ISLocalToGlobalMapping ltogi; /* Map from global to local indices of the i-th outer local subdomain. */ 4266a4f0f73SDmitry Karpeev PetscInt *ioidxi; /* Local indices of the i-th local inner subdomain within the local outer subdomain. */ 4276a4f0f73SDmitry Karpeev PetscInt ioni; /* Number of indices in ioidxi; if ioni != ini the inner subdomain is not a subdomain of the outer subdomain (error). */ 4286a4f0f73SDmitry Karpeev PetscInt k; 4296a4f0f73SDmitry Karpeev ierr = ISGetLocalSize(osm->iis[i],&ini);CHKERRQ(ierr); 4306a4f0f73SDmitry Karpeev ierr = ISGetLocalSize(osm->ois[i],&oni);CHKERRQ(ierr); 4316a4f0f73SDmitry Karpeev ierr = ISGetIndices(osm->iis[i],&iidxi);CHKERRQ(ierr); 4326a4f0f73SDmitry Karpeev ierr = PetscMemcpy(iidx+in, iidxi, sizeof(PetscInt)*ini);CHKERRQ(ierr); 4336a4f0f73SDmitry Karpeev ierr = ISLocalToGlobalMappingCreateIS(osm->ois[i],<ogi);CHKERRQ(ierr); 4346a4f0f73SDmitry Karpeev ioidxi = ioidx+in; 4356a4f0f73SDmitry Karpeev ierr = ISGlobalToLocalMappingApply(ltogi,IS_GTOLM_DROP,ini,iidxi,&ioni,ioidxi);CHKERRQ(ierr); 4366a4f0f73SDmitry Karpeev ierr = ISLocalToGlobalMappingDestroy(<ogi);CHKERRQ(ierr); 4376a4f0f73SDmitry Karpeev ierr = ISRestoreIndices(osm->iis[i], &iidxi);CHKERRQ(ierr); 4386a4f0f73SDmitry Karpeev if (ioni != ini) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Inner subdomain %D contains %D indices outside of its outer subdomain", i, ini - ioni); 4392fa5cd67SKarl Rupp for (k = 0; k < ini; ++k) ioidxi[k] += gofirst+on; 4406a4f0f73SDmitry Karpeev in += ini; 4416a4f0f73SDmitry Karpeev on += oni; 442f746d493SDmitry Karpeev } 443ce94432eSBarry Smith ierr = ISCreateGeneral(PetscObjectComm((PetscObject)pc), in, iidx, PETSC_OWN_POINTER, &giis);CHKERRQ(ierr); 444ce94432eSBarry Smith ierr = ISCreateGeneral(PetscObjectComm((PetscObject)pc), in, ioidx, PETSC_OWN_POINTER, &giois);CHKERRQ(ierr); 4456a4f0f73SDmitry Karpeev ierr = VecScatterCreate(y,giis,osm->gy,giois,&osm->girestriction);CHKERRQ(ierr); 4466a4f0f73SDmitry Karpeev ierr = VecDestroy(&y);CHKERRQ(ierr); 4476a4f0f73SDmitry Karpeev ierr = ISDestroy(&giis);CHKERRQ(ierr); 4486a4f0f73SDmitry Karpeev ierr = ISDestroy(&giois);CHKERRQ(ierr); 449b1a0a8a3SJed Brown } 4506a4f0f73SDmitry Karpeev ierr = ISDestroy(&goid);CHKERRQ(ierr); 4512fa5cd67SKarl Rupp 452f746d493SDmitry Karpeev /* Create the subdomain work vectors. */ 453785e854fSJed Brown ierr = PetscMalloc1(osm->n,&osm->x);CHKERRQ(ierr); 454785e854fSJed Brown ierr = PetscMalloc1(osm->n,&osm->y);CHKERRQ(ierr); 455f746d493SDmitry Karpeev ierr = VecGetArray(osm->gx, &gxarray);CHKERRQ(ierr); 456f746d493SDmitry Karpeev ierr = VecGetArray(osm->gy, &gyarray);CHKERRQ(ierr); 4576a4f0f73SDmitry Karpeev for (i=0, on=0; i<osm->n; ++i, on += oni) { 4586a4f0f73SDmitry Karpeev PetscInt oNi; 4596a4f0f73SDmitry Karpeev ierr = ISGetLocalSize(osm->ois[i],&oni);CHKERRQ(ierr); 4606a4f0f73SDmitry Karpeev ierr = ISGetSize(osm->ois[i],&oNi);CHKERRQ(ierr); 4616a4f0f73SDmitry Karpeev ierr = VecCreateMPIWithArray(((PetscObject)(osm->ois[i]))->comm,1,oni,oNi,gxarray+on,&osm->x[i]);CHKERRQ(ierr); 4626a4f0f73SDmitry Karpeev ierr = VecCreateMPIWithArray(((PetscObject)(osm->ois[i]))->comm,1,oni,oNi,gyarray+on,&osm->y[i]);CHKERRQ(ierr); 463b1a0a8a3SJed Brown } 464f746d493SDmitry Karpeev ierr = VecRestoreArray(osm->gx, &gxarray);CHKERRQ(ierr); 465f746d493SDmitry Karpeev ierr = VecRestoreArray(osm->gy, &gyarray);CHKERRQ(ierr); 466*8f3b4b4dSDmitry Karpeev /* Create the subdomain solvers */ 467785e854fSJed Brown ierr = PetscMalloc1(osm->n,&osm->ksp);CHKERRQ(ierr); 468*8f3b4b4dSDmitry Karpeev for (i=0; i<osm->n; i++) { 469*8f3b4b4dSDmitry Karpeev char subprefix[PETSC_MAX_PATH_LEN+1]; 4706a4f0f73SDmitry Karpeev ierr = KSPCreate(((PetscObject)(osm->ois[i]))->comm,&ksp);CHKERRQ(ierr); 4713bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ksp);CHKERRQ(ierr); 472b1a0a8a3SJed Brown ierr = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,1);CHKERRQ(ierr); 473b1a0a8a3SJed Brown ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr); 474*8f3b4b4dSDmitry Karpeev ierr = KSPGetPC(ksp,&subpc);CHKERRQ(ierr); /* Why do we need this here? */ 475*8f3b4b4dSDmitry Karpeev if (subdomain_dm) { 476*8f3b4b4dSDmitry Karpeev ierr = KSPSetDM(ksp,subdomain_dm[i]);CHKERRQ(ierr); 477*8f3b4b4dSDmitry Karpeev ierr = DMDestroy(subdomain_dm+i);CHKERRQ(ierr); 478*8f3b4b4dSDmitry Karpeev } 479b1a0a8a3SJed Brown ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 480b1a0a8a3SJed Brown ierr = KSPSetOptionsPrefix(ksp,prefix);CHKERRQ(ierr); 481*8f3b4b4dSDmitry Karpeev if (subdomain_names && subdomain_names[i]) { 482*8f3b4b4dSDmitry Karpeev ierr = PetscSNPrintf(subprefix,PETSC_MAX_PATH_LEN,"sub_%s_",subdomain_names[i]);CHKERRQ(ierr); 483*8f3b4b4dSDmitry Karpeev ierr = KSPAppendOptionsPrefix(ksp,subprefix);CHKERRQ(ierr); 484*8f3b4b4dSDmitry Karpeev ierr = PetscFree(subdomain_names[i]);CHKERRQ(ierr); 485*8f3b4b4dSDmitry Karpeev } 486b1a0a8a3SJed Brown ierr = KSPAppendOptionsPrefix(ksp,"sub_");CHKERRQ(ierr); 487b1a0a8a3SJed Brown osm->ksp[i] = ksp; 488b1a0a8a3SJed Brown } 489*8f3b4b4dSDmitry Karpeev ierr = PetscFree(subdomain_dm);CHKERRQ(ierr); 490*8f3b4b4dSDmitry Karpeev ierr = PetscFree(subdomain_names);CHKERRQ(ierr); 491b1a0a8a3SJed Brown scall = MAT_INITIAL_MATRIX; 492b1a0a8a3SJed Brown 493*8f3b4b4dSDmitry Karpeev } else { /*if (pc->setupcalled)*/ 494b1a0a8a3SJed Brown /* 495*8f3b4b4dSDmitry Karpeev Destroy the submatrices from the previous iteration 496b1a0a8a3SJed Brown */ 497b1a0a8a3SJed Brown if (pc->flag == DIFFERENT_NONZERO_PATTERN) { 498f746d493SDmitry Karpeev ierr = MatDestroyMatrices(osm->n,&osm->pmat);CHKERRQ(ierr); 499b1a0a8a3SJed Brown scall = MAT_INITIAL_MATRIX; 500b1a0a8a3SJed Brown } 501b1a0a8a3SJed Brown } 502b1a0a8a3SJed Brown 503b1a0a8a3SJed Brown /* 504f746d493SDmitry Karpeev Extract out the submatrices. 505b1a0a8a3SJed Brown */ 50681a5c4aaSDmitry Karpeev if (size > 1) { 507*8f3b4b4dSDmitry Karpeev ierr = MatGetSubMatricesMPI(pc->pmat,osm->n,osm->ois,osm->ois,scall,&osm->pmat);CHKERRQ(ierr); 5082fa5cd67SKarl Rupp } else { 5096a4f0f73SDmitry Karpeev ierr = MatGetSubMatrices(pc->pmat,osm->n,osm->ois,osm->ois,scall,&osm->pmat);CHKERRQ(ierr); 51081a5c4aaSDmitry Karpeev } 511b1a0a8a3SJed Brown if (scall == MAT_INITIAL_MATRIX) { 512b1a0a8a3SJed Brown ierr = PetscObjectGetOptionsPrefix((PetscObject)pc->pmat,&pprefix);CHKERRQ(ierr); 513f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 5143bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)osm->pmat[i]);CHKERRQ(ierr); 515b1a0a8a3SJed Brown ierr = PetscObjectSetOptionsPrefix((PetscObject)osm->pmat[i],pprefix);CHKERRQ(ierr); 516b1a0a8a3SJed Brown } 517b1a0a8a3SJed Brown } 518b1a0a8a3SJed Brown 519b1a0a8a3SJed Brown /* Return control to the user so that the submatrices can be modified (e.g., to apply 520b1a0a8a3SJed Brown different boundary conditions for the submatrices than for the global problem) */ 5216a4f0f73SDmitry Karpeev ierr = PCModifySubMatrices(pc,osm->n,osm->ois,osm->ois,osm->pmat,pc->modifysubmatricesP);CHKERRQ(ierr); 522b1a0a8a3SJed Brown 523b1a0a8a3SJed Brown /* 5246a4f0f73SDmitry Karpeev Loop over submatrices putting them into local ksps 525b1a0a8a3SJed Brown */ 526f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 52723ee1639SBarry Smith ierr = KSPSetOperators(osm->ksp[i],osm->pmat[i],osm->pmat[i]);CHKERRQ(ierr); 528b1a0a8a3SJed Brown if (!pc->setupcalled) { 529b1a0a8a3SJed Brown ierr = KSPSetFromOptions(osm->ksp[i]);CHKERRQ(ierr); 530b1a0a8a3SJed Brown } 531b1a0a8a3SJed Brown } 532b1a0a8a3SJed Brown PetscFunctionReturn(0); 533b1a0a8a3SJed Brown } 534b1a0a8a3SJed Brown 535b1a0a8a3SJed Brown #undef __FUNCT__ 536f746d493SDmitry Karpeev #define __FUNCT__ "PCSetUpOnBlocks_GASM" 537f746d493SDmitry Karpeev static PetscErrorCode PCSetUpOnBlocks_GASM(PC pc) 538b1a0a8a3SJed Brown { 539f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 540b1a0a8a3SJed Brown PetscErrorCode ierr; 541b1a0a8a3SJed Brown PetscInt i; 542b1a0a8a3SJed Brown 543b1a0a8a3SJed Brown PetscFunctionBegin; 544f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 545b1a0a8a3SJed Brown ierr = KSPSetUp(osm->ksp[i]);CHKERRQ(ierr); 546b1a0a8a3SJed Brown } 547b1a0a8a3SJed Brown PetscFunctionReturn(0); 548b1a0a8a3SJed Brown } 549b1a0a8a3SJed Brown 550b1a0a8a3SJed Brown #undef __FUNCT__ 551f746d493SDmitry Karpeev #define __FUNCT__ "PCApply_GASM" 552f746d493SDmitry Karpeev static PetscErrorCode PCApply_GASM(PC pc,Vec x,Vec y) 553b1a0a8a3SJed Brown { 554f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 555b1a0a8a3SJed Brown PetscErrorCode ierr; 556f746d493SDmitry Karpeev PetscInt i; 557b1a0a8a3SJed Brown ScatterMode forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE; 558b1a0a8a3SJed Brown 559b1a0a8a3SJed Brown PetscFunctionBegin; 560b1a0a8a3SJed Brown /* 5616a4f0f73SDmitry Karpeev Support for limiting the restriction or interpolation only to the inner 562b1a0a8a3SJed Brown subdomain values (leaving the other values 0). 563b1a0a8a3SJed Brown */ 564f746d493SDmitry Karpeev if (!(osm->type & PC_GASM_RESTRICT)) { 565b1a0a8a3SJed Brown /* have to zero the work RHS since scatter may leave some slots empty */ 566f746d493SDmitry Karpeev ierr = VecZeroEntries(osm->gx);CHKERRQ(ierr); 5676a4f0f73SDmitry Karpeev ierr = VecScatterBegin(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 5682fa5cd67SKarl Rupp } else { 5696a4f0f73SDmitry Karpeev ierr = VecScatterBegin(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 570b1a0a8a3SJed Brown } 5716a4f0f73SDmitry Karpeev ierr = VecZeroEntries(osm->gy);CHKERRQ(ierr); 5726a4f0f73SDmitry Karpeev if (!(osm->type & PC_GASM_RESTRICT)) { 5736a4f0f73SDmitry Karpeev ierr = VecScatterEnd(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 5742fa5cd67SKarl Rupp } else { 5756a4f0f73SDmitry Karpeev ierr = VecScatterEnd(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 5766a4f0f73SDmitry Karpeev } 57786b96d36SDmitry Karpeev /* do the subdomain solves */ 57886b96d36SDmitry Karpeev for (i=0; i<osm->n; ++i) { 579b1a0a8a3SJed Brown ierr = KSPSolve(osm->ksp[i],osm->x[i],osm->y[i]);CHKERRQ(ierr); 580b1a0a8a3SJed Brown } 5816a4f0f73SDmitry Karpeev ierr = VecZeroEntries(y);CHKERRQ(ierr); 5826a4f0f73SDmitry Karpeev if (!(osm->type & PC_GASM_INTERPOLATE)) { 5836a4f0f73SDmitry Karpeev ierr = VecScatterBegin(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); 5846a4f0f73SDmitry Karpeev ierr = VecScatterEnd(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); PetscFunctionReturn(0); 5852fa5cd67SKarl Rupp } else { 5866a4f0f73SDmitry Karpeev ierr = VecScatterBegin(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); 5876a4f0f73SDmitry Karpeev ierr = VecScatterEnd(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); PetscFunctionReturn(0); 5886a4f0f73SDmitry Karpeev } 589b1a0a8a3SJed Brown } 590b1a0a8a3SJed Brown 591b1a0a8a3SJed Brown #undef __FUNCT__ 592f746d493SDmitry Karpeev #define __FUNCT__ "PCApplyTranspose_GASM" 593f746d493SDmitry Karpeev static PetscErrorCode PCApplyTranspose_GASM(PC pc,Vec x,Vec y) 594b1a0a8a3SJed Brown { 595f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 596b1a0a8a3SJed Brown PetscErrorCode ierr; 597f746d493SDmitry Karpeev PetscInt i; 598b1a0a8a3SJed Brown ScatterMode forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE; 599b1a0a8a3SJed Brown 600b1a0a8a3SJed Brown PetscFunctionBegin; 601b1a0a8a3SJed Brown /* 602b1a0a8a3SJed Brown Support for limiting the restriction or interpolation to only local 603b1a0a8a3SJed Brown subdomain values (leaving the other values 0). 604b1a0a8a3SJed Brown 605f746d493SDmitry Karpeev Note: these are reversed from the PCApply_GASM() because we are applying the 606b1a0a8a3SJed Brown transpose of the three terms 607b1a0a8a3SJed Brown */ 608f746d493SDmitry Karpeev if (!(osm->type & PC_GASM_INTERPOLATE)) { 609b1a0a8a3SJed Brown /* have to zero the work RHS since scatter may leave some slots empty */ 610f746d493SDmitry Karpeev ierr = VecZeroEntries(osm->gx);CHKERRQ(ierr); 6116a4f0f73SDmitry Karpeev ierr = VecScatterBegin(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 6122fa5cd67SKarl Rupp } else { 6136a4f0f73SDmitry Karpeev ierr = VecScatterBegin(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 614b1a0a8a3SJed Brown } 6156a4f0f73SDmitry Karpeev ierr = VecZeroEntries(osm->gy);CHKERRQ(ierr); 6166a4f0f73SDmitry Karpeev if (!(osm->type & PC_GASM_INTERPOLATE)) { 6176a4f0f73SDmitry Karpeev ierr = VecScatterEnd(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 6182fa5cd67SKarl Rupp } else { 6196a4f0f73SDmitry Karpeev ierr = VecScatterEnd(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 6206a4f0f73SDmitry Karpeev } 621b1a0a8a3SJed Brown /* do the local solves */ 622ab3e923fSDmitry 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. */ 623b1a0a8a3SJed Brown ierr = KSPSolveTranspose(osm->ksp[i],osm->x[i],osm->y[i]);CHKERRQ(ierr); 624b1a0a8a3SJed Brown } 6256a4f0f73SDmitry Karpeev ierr = VecZeroEntries(y);CHKERRQ(ierr); 6266a4f0f73SDmitry Karpeev if (!(osm->type & PC_GASM_RESTRICT)) { 6276a4f0f73SDmitry Karpeev ierr = VecScatterBegin(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); 6286a4f0f73SDmitry Karpeev ierr = VecScatterEnd(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); 6292fa5cd67SKarl Rupp } else { 6306a4f0f73SDmitry Karpeev ierr = VecScatterBegin(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); 6316a4f0f73SDmitry Karpeev ierr = VecScatterEnd(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); 6326a4f0f73SDmitry Karpeev } 633b1a0a8a3SJed Brown PetscFunctionReturn(0); 634b1a0a8a3SJed Brown } 635b1a0a8a3SJed Brown 636b1a0a8a3SJed Brown #undef __FUNCT__ 637a06653b4SBarry Smith #define __FUNCT__ "PCReset_GASM" 638a06653b4SBarry Smith static PetscErrorCode PCReset_GASM(PC pc) 639b1a0a8a3SJed Brown { 640f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 641b1a0a8a3SJed Brown PetscErrorCode ierr; 642b1a0a8a3SJed Brown PetscInt i; 643b1a0a8a3SJed Brown 644b1a0a8a3SJed Brown PetscFunctionBegin; 645b1a0a8a3SJed Brown if (osm->ksp) { 646f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 647a06653b4SBarry Smith ierr = KSPReset(osm->ksp[i]);CHKERRQ(ierr); 648b1a0a8a3SJed Brown } 649b1a0a8a3SJed Brown } 650b1a0a8a3SJed Brown if (osm->pmat) { 651f746d493SDmitry Karpeev if (osm->n > 0) { 652ab3e923fSDmitry Karpeev ierr = MatDestroyMatrices(osm->n,&osm->pmat);CHKERRQ(ierr); 653b1a0a8a3SJed Brown } 654b1a0a8a3SJed Brown } 655d34fcf5fSBarry Smith if (osm->x) { 656f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 657fcfd50ebSBarry Smith ierr = VecDestroy(&osm->x[i]);CHKERRQ(ierr); 658fcfd50ebSBarry Smith ierr = VecDestroy(&osm->y[i]);CHKERRQ(ierr); 659b1a0a8a3SJed Brown } 660d34fcf5fSBarry Smith } 6616bf464f9SBarry Smith ierr = VecDestroy(&osm->gx);CHKERRQ(ierr); 6626bf464f9SBarry Smith ierr = VecDestroy(&osm->gy);CHKERRQ(ierr); 663ab3e923fSDmitry Karpeev 6646a4f0f73SDmitry Karpeev ierr = VecScatterDestroy(&osm->gorestriction);CHKERRQ(ierr); 6656a4f0f73SDmitry Karpeev ierr = VecScatterDestroy(&osm->girestriction);CHKERRQ(ierr); 666*8f3b4b4dSDmitry Karpeev if (!osm->user_subdomains) { 6676a4f0f73SDmitry Karpeev ierr = PCGASMDestroySubdomains(osm->n,osm->ois,osm->iis);CHKERRQ(ierr); 668*8f3b4b4dSDmitry Karpeev osm->N = PETSC_DETERMINE; 669*8f3b4b4dSDmitry Karpeev osm->nmax = PETSC_DETERMINE; 6706a4f0f73SDmitry Karpeev osm->ois = 0; 6716a4f0f73SDmitry Karpeev osm->iis = 0; 672*8f3b4b4dSDmitry Karpeev } 673a06653b4SBarry Smith PetscFunctionReturn(0); 674a06653b4SBarry Smith } 675a06653b4SBarry Smith 676a06653b4SBarry Smith #undef __FUNCT__ 677a06653b4SBarry Smith #define __FUNCT__ "PCDestroy_GASM" 678a06653b4SBarry Smith static PetscErrorCode PCDestroy_GASM(PC pc) 679a06653b4SBarry Smith { 680a06653b4SBarry Smith PC_GASM *osm = (PC_GASM*)pc->data; 681a06653b4SBarry Smith PetscErrorCode ierr; 682a06653b4SBarry Smith PetscInt i; 683a06653b4SBarry Smith 684a06653b4SBarry Smith PetscFunctionBegin; 685a06653b4SBarry Smith ierr = PCReset_GASM(pc);CHKERRQ(ierr); 686a06653b4SBarry Smith if (osm->ksp) { 687a06653b4SBarry Smith for (i=0; i<osm->n; i++) { 6886bf464f9SBarry Smith ierr = KSPDestroy(&osm->ksp[i]);CHKERRQ(ierr); 689a06653b4SBarry Smith } 690a06653b4SBarry Smith ierr = PetscFree(osm->ksp);CHKERRQ(ierr); 691a06653b4SBarry Smith } 692a06653b4SBarry Smith ierr = PetscFree(osm->x);CHKERRQ(ierr); 693a06653b4SBarry Smith ierr = PetscFree(osm->y);CHKERRQ(ierr); 694c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 695b1a0a8a3SJed Brown PetscFunctionReturn(0); 696b1a0a8a3SJed Brown } 697b1a0a8a3SJed Brown 698b1a0a8a3SJed Brown #undef __FUNCT__ 699f746d493SDmitry Karpeev #define __FUNCT__ "PCSetFromOptions_GASM" 7008c34d3f5SBarry Smith static PetscErrorCode PCSetFromOptions_GASM(PetscOptions *PetscOptionsObject,PC pc) 701a6dfd86eSKarl Rupp { 702f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 703b1a0a8a3SJed Brown PetscErrorCode ierr; 704b1a0a8a3SJed Brown PetscInt blocks,ovl; 705b1a0a8a3SJed Brown PetscBool symset,flg; 706f746d493SDmitry Karpeev PCGASMType gasmtype; 707b1a0a8a3SJed Brown 708b1a0a8a3SJed Brown PetscFunctionBegin; 709b1a0a8a3SJed Brown /* set the type to symmetric if matrix is symmetric */ 710b1a0a8a3SJed Brown if (!osm->type_set && pc->pmat) { 711b1a0a8a3SJed Brown ierr = MatIsSymmetricKnown(pc->pmat,&symset,&flg);CHKERRQ(ierr); 7122fa5cd67SKarl Rupp if (symset && flg) osm->type = PC_GASM_BASIC; 713b1a0a8a3SJed Brown } 714e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"Generalized additive Schwarz options");CHKERRQ(ierr); 715*8f3b4b4dSDmitry 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); 716*8f3b4b4dSDmitry Karpeev ierr = PetscOptionsInt("-pc_gasm_total_subdomains","Total number of subdomains across communicator","PCGASMSetTotalSubdomains",osm->N,&blocks,&flg);CHKERRQ(ierr); 71765db9045SDmitry Karpeev if (flg) { 718*8f3b4b4dSDmitry Karpeev ierr = PCGASMSetTotalSubdomains(pc,blocks);CHKERRQ(ierr); 71965db9045SDmitry Karpeev } 72006b43e7eSDmitry Karpeev ierr = PetscOptionsInt("-pc_gasm_overlap","Number of overlapping degrees of freedom","PCGASMSetOverlap",osm->overlap,&ovl,&flg);CHKERRQ(ierr); 72165db9045SDmitry Karpeev if (flg) { 72265db9045SDmitry Karpeev ierr = PCGASMSetOverlap(pc,ovl);CHKERRQ(ierr); 723d709ea83SDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 72465db9045SDmitry Karpeev } 725b1a0a8a3SJed Brown flg = PETSC_FALSE; 726f746d493SDmitry Karpeev ierr = PetscOptionsEnum("-pc_gasm_type","Type of restriction/extension","PCGASMSetType",PCGASMTypes,(PetscEnum)osm->type,(PetscEnum*)&gasmtype,&flg);CHKERRQ(ierr); 727f746d493SDmitry Karpeev if (flg) {ierr = PCGASMSetType(pc,gasmtype);CHKERRQ(ierr);} 728b1a0a8a3SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 729b1a0a8a3SJed Brown PetscFunctionReturn(0); 730b1a0a8a3SJed Brown } 731b1a0a8a3SJed Brown 732b1a0a8a3SJed Brown /*------------------------------------------------------------------------------------*/ 733b1a0a8a3SJed Brown 734b1a0a8a3SJed Brown #undef __FUNCT__ 735*8f3b4b4dSDmitry Karpeev #define __FUNCT__ "PCGASMSetTotalSubdomains" 736*8f3b4b4dSDmitry Karpeev /*@ 737*8f3b4b4dSDmitry Karpeev PCGASMSetTotalSubdomains - sets the total number of subdomains to use across the 738*8f3b4b4dSDmitry Karpeev communicator. 739*8f3b4b4dSDmitry Karpeev Logically collective on pc 740*8f3b4b4dSDmitry Karpeev 741*8f3b4b4dSDmitry Karpeev Input Parameters: 742*8f3b4b4dSDmitry Karpeev + pc - the preconditioner 743*8f3b4b4dSDmitry Karpeev - N - total number of subdomains 744*8f3b4b4dSDmitry Karpeev 745*8f3b4b4dSDmitry Karpeev 746*8f3b4b4dSDmitry Karpeev Level: beginner 747*8f3b4b4dSDmitry Karpeev 748*8f3b4b4dSDmitry Karpeev .keywords: PC, ASM, DM, set, subdomains, additive Schwarz 749*8f3b4b4dSDmitry Karpeev 750*8f3b4b4dSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMSetOverlap() 751*8f3b4b4dSDmitry Karpeev PCGASMCreateSubdomains2D() 752*8f3b4b4dSDmitry Karpeev @*/ 753*8f3b4b4dSDmitry Karpeev PetscErrorCode PCGASMSetTotalSubdomains(PC pc,PetscInt N) 754*8f3b4b4dSDmitry Karpeev { 755*8f3b4b4dSDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 756*8f3b4b4dSDmitry Karpeev PetscMPIInt size,rank; 757*8f3b4b4dSDmitry Karpeev PetscErrorCode ierr; 758*8f3b4b4dSDmitry Karpeev 759*8f3b4b4dSDmitry Karpeev PetscFunctionBegin; 760*8f3b4b4dSDmitry 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); 761*8f3b4b4dSDmitry Karpeev if (pc->setupcalled) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetTotalSubdomains() should be called before calling PCSetUp()."); 762*8f3b4b4dSDmitry Karpeev 763*8f3b4b4dSDmitry Karpeev ierr = PCGASMDestroySubdomains(osm->n,osm->iis,osm->ois);CHKERRQ(ierr); 764*8f3b4b4dSDmitry Karpeev osm->ois = osm->iis = NULL; 765*8f3b4b4dSDmitry Karpeev 766*8f3b4b4dSDmitry Karpeev ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr); 767*8f3b4b4dSDmitry Karpeev ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);CHKERRQ(ierr); 768*8f3b4b4dSDmitry Karpeev osm->N = N; 769*8f3b4b4dSDmitry Karpeev osm->n = PETSC_DETERMINE; 770*8f3b4b4dSDmitry Karpeev osm->nmax = PETSC_DETERMINE; 771*8f3b4b4dSDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 772*8f3b4b4dSDmitry Karpeev PetscFunctionReturn(0); 773*8f3b4b4dSDmitry Karpeev } 774*8f3b4b4dSDmitry Karpeev 775*8f3b4b4dSDmitry Karpeev #undef __FUNCT__ 77606b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMSetSubdomains_GASM" 777*8f3b4b4dSDmitry Karpeev static PetscErrorCode PCGASMSetSubdomains_GASM(PC pc,PetscMPIInt n,IS iis[],IS ois[]) 778b1a0a8a3SJed Brown { 779f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 780b1a0a8a3SJed Brown PetscErrorCode ierr; 781b1a0a8a3SJed Brown PetscInt i; 782b1a0a8a3SJed Brown 783b1a0a8a3SJed Brown PetscFunctionBegin; 784*8f3b4b4dSDmitry Karpeev if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Each process must have 1 or more subdomains, got n = %D",n); 785*8f3b4b4dSDmitry Karpeev if (pc->setupcalled) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetSubdomains() should be called before calling PCSetUp()."); 786b1a0a8a3SJed Brown 7876a4f0f73SDmitry Karpeev ierr = PCGASMDestroySubdomains(osm->n,osm->iis,osm->ois);CHKERRQ(ierr); 788*8f3b4b4dSDmitry Karpeev osm->iis = osm->ois = NULL; 789*8f3b4b4dSDmitry Karpeev osm->n = n; 790*8f3b4b4dSDmitry Karpeev osm->N = PETSC_DETERMINE; 791*8f3b4b4dSDmitry Karpeev osm->nmax = PETSC_DETERMINE; 792a35b7b57SDmitry Karpeev if (ois) { 793785e854fSJed Brown ierr = PetscMalloc1(n,&osm->ois);CHKERRQ(ierr); 794*8f3b4b4dSDmitry Karpeev for (i=0; i<n; i++) { 795*8f3b4b4dSDmitry Karpeev ierr = PetscObjectReference((PetscObject)ois[i]);CHKERRQ(ierr); 796*8f3b4b4dSDmitry Karpeev osm->ois[i] = ois[i]; 797*8f3b4b4dSDmitry Karpeev } 798*8f3b4b4dSDmitry Karpeev /* 799*8f3b4b4dSDmitry Karpeev Since the user set the outer subdomains, even if nontrivial overlap was requested via PCGASMSetOverlap(), 800*8f3b4b4dSDmitry Karpeev it will be ignored. To avoid confusion later on (e.g., when viewing the PC), the overlap size is set to -1. 801*8f3b4b4dSDmitry Karpeev */ 802b1a0a8a3SJed Brown osm->overlap = -1; 803a35b7b57SDmitry Karpeev if (!iis) { 804785e854fSJed Brown ierr = PetscMalloc1(n,&osm->iis);CHKERRQ(ierr); 805a35b7b57SDmitry Karpeev for (i=0; i<n; i++) { 806*8f3b4b4dSDmitry Karpeev for (i=0; i<n; i++) { 8070e8f3646SDmitry Karpeev ierr = PetscObjectReference((PetscObject)ois[i]);CHKERRQ(ierr); 808a35b7b57SDmitry Karpeev osm->iis[i] = ois[i]; 809a35b7b57SDmitry Karpeev } 810a35b7b57SDmitry Karpeev } 811a35b7b57SDmitry Karpeev } 812*8f3b4b4dSDmitry Karpeev } 813a35b7b57SDmitry Karpeev if (iis) { 814785e854fSJed Brown ierr = PetscMalloc1(n,&osm->iis);CHKERRQ(ierr); 815*8f3b4b4dSDmitry Karpeev for (i=0; i<n; i++) { 816*8f3b4b4dSDmitry Karpeev ierr = PetscObjectReference((PetscObject)iis[i]);CHKERRQ(ierr); 817*8f3b4b4dSDmitry Karpeev osm->iis[i] = iis[i]; 818*8f3b4b4dSDmitry Karpeev } 819a35b7b57SDmitry Karpeev if (!ois) { 820785e854fSJed Brown ierr = PetscMalloc1(n,&osm->ois);CHKERRQ(ierr); 821a35b7b57SDmitry Karpeev for (i=0; i<n; i++) { 822*8f3b4b4dSDmitry Karpeev for (i=0; i<n; i++) { 823a35b7b57SDmitry Karpeev ierr = PetscObjectReference((PetscObject)iis[i]);CHKERRQ(ierr); 824a35b7b57SDmitry Karpeev osm->ois[i] = iis[i]; 825a35b7b57SDmitry Karpeev } 826*8f3b4b4dSDmitry Karpeev } 827*8f3b4b4dSDmitry Karpeev /* If necessary, osm->ois[i] will be expanded using the requested overlap size in PCSetUp_GASM(). */ 828a35b7b57SDmitry Karpeev } 829a35b7b57SDmitry Karpeev } 830*8f3b4b4dSDmitry Karpeev if (iis) osm->user_subdomains = PETSC_TRUE; 831b1a0a8a3SJed Brown PetscFunctionReturn(0); 832b1a0a8a3SJed Brown } 833b1a0a8a3SJed Brown 834b1a0a8a3SJed Brown 835b1a0a8a3SJed Brown #undef __FUNCT__ 836f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetOverlap_GASM" 837f7a08781SBarry Smith static PetscErrorCode PCGASMSetOverlap_GASM(PC pc,PetscInt ovl) 838b1a0a8a3SJed Brown { 839f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 840b1a0a8a3SJed Brown 841b1a0a8a3SJed Brown PetscFunctionBegin; 842ce94432eSBarry Smith if (ovl < 0) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap value requested"); 843ce94432eSBarry Smith if (pc->setupcalled && ovl != osm->overlap) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetOverlap() should be called before PCSetUp()."); 8442fa5cd67SKarl Rupp if (!pc->setupcalled) osm->overlap = ovl; 845b1a0a8a3SJed Brown PetscFunctionReturn(0); 846b1a0a8a3SJed Brown } 847b1a0a8a3SJed Brown 848b1a0a8a3SJed Brown #undef __FUNCT__ 849f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetType_GASM" 850f7a08781SBarry Smith static PetscErrorCode PCGASMSetType_GASM(PC pc,PCGASMType type) 851b1a0a8a3SJed Brown { 852f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 853b1a0a8a3SJed Brown 854b1a0a8a3SJed Brown PetscFunctionBegin; 855b1a0a8a3SJed Brown osm->type = type; 856b1a0a8a3SJed Brown osm->type_set = PETSC_TRUE; 857b1a0a8a3SJed Brown PetscFunctionReturn(0); 858b1a0a8a3SJed Brown } 859b1a0a8a3SJed Brown 860b1a0a8a3SJed Brown #undef __FUNCT__ 861f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetSortIndices_GASM" 862f7a08781SBarry Smith static PetscErrorCode PCGASMSetSortIndices_GASM(PC pc,PetscBool doSort) 863b1a0a8a3SJed Brown { 864f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 865b1a0a8a3SJed Brown 866b1a0a8a3SJed Brown PetscFunctionBegin; 867b1a0a8a3SJed Brown osm->sort_indices = doSort; 868b1a0a8a3SJed Brown PetscFunctionReturn(0); 869b1a0a8a3SJed Brown } 870b1a0a8a3SJed Brown 871b1a0a8a3SJed Brown #undef __FUNCT__ 872f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMGetSubKSP_GASM" 87344fe18e5SDmitry Karpeev /* 874*8f3b4b4dSDmitry Karpeev FIXME: This routine might need to be modified now that multiple ranks per subdomain are allowed. 87544fe18e5SDmitry Karpeev In particular, it would upset the global subdomain number calculation. 87644fe18e5SDmitry Karpeev */ 877f7a08781SBarry Smith static PetscErrorCode PCGASMGetSubKSP_GASM(PC pc,PetscInt *n,PetscInt *first,KSP **ksp) 878b1a0a8a3SJed Brown { 879f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 880b1a0a8a3SJed Brown PetscErrorCode ierr; 881b1a0a8a3SJed Brown 882b1a0a8a3SJed Brown PetscFunctionBegin; 883ce94432eSBarry 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"); 884b1a0a8a3SJed Brown 8852fa5cd67SKarl Rupp if (n) *n = osm->n; 886ab3e923fSDmitry Karpeev if (first) { 887ce94432eSBarry Smith ierr = MPI_Scan(&osm->n,first,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr); 888ab3e923fSDmitry Karpeev *first -= osm->n; 889b1a0a8a3SJed Brown } 890b1a0a8a3SJed Brown if (ksp) { 891b1a0a8a3SJed Brown /* Assume that local solves are now different; not necessarily 89206b43e7eSDmitry Karpeev true, though! This flag is used only for PCView_GASM() */ 893b1a0a8a3SJed Brown *ksp = osm->ksp; 8946a4f0f73SDmitry Karpeev osm->same_subdomain_solvers = PETSC_FALSE; 895b1a0a8a3SJed Brown } 896b1a0a8a3SJed Brown PetscFunctionReturn(0); 897ab3e923fSDmitry Karpeev } /* PCGASMGetSubKSP_GASM() */ 898b1a0a8a3SJed Brown 899b1a0a8a3SJed Brown #undef __FUNCT__ 90006b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMSetSubdomains" 901b1a0a8a3SJed Brown /*@C 90206b43e7eSDmitry Karpeev PCGASMSetSubdomains - Sets the subdomains for this processor 90306b43e7eSDmitry Karpeev for the additive Schwarz preconditioner. 904b1a0a8a3SJed Brown 905*8f3b4b4dSDmitry Karpeev Collective on pc 906b1a0a8a3SJed Brown 907b1a0a8a3SJed Brown Input Parameters: 908*8f3b4b4dSDmitry Karpeev + pc - the preconditioner object 90906b43e7eSDmitry Karpeev . n - the number of subdomains for this processor 910*8f3b4b4dSDmitry Karpeev . iis - the index sets that define the inner subdomains (or NULL for PETSc to determine subdomains) 911*8f3b4b4dSDmitry 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) 912b1a0a8a3SJed Brown 913b1a0a8a3SJed Brown Notes: 91406b43e7eSDmitry Karpeev The IS indices use the parallel, global numbering of the vector entries. 9156a4f0f73SDmitry Karpeev Inner subdomains are those where the correction is applied. 9166a4f0f73SDmitry Karpeev Outer subdomains are those where the residual necessary to obtain the 9176a4f0f73SDmitry Karpeev corrections is obtained (see PCGASMType for the use of inner/outer subdomains). 9186a4f0f73SDmitry Karpeev Both inner and outer subdomains can extend over several processors. 9196a4f0f73SDmitry Karpeev This processor's portion of a subdomain is known as a local subdomain. 9206a4f0f73SDmitry Karpeev 9216a4f0f73SDmitry Karpeev By default the GASM preconditioner uses 1 (local) subdomain per processor. 9226a4f0f73SDmitry Karpeev 923b1a0a8a3SJed Brown 924b1a0a8a3SJed Brown Level: advanced 925b1a0a8a3SJed Brown 92606b43e7eSDmitry Karpeev .keywords: PC, GASM, set, subdomains, additive Schwarz 927b1a0a8a3SJed Brown 928*8f3b4b4dSDmitry Karpeev .seealso: PCGASMSetNumSubdomains(), PCGASMSetOverlap(), PCGASMGetSubKSP(), 92906b43e7eSDmitry Karpeev PCGASMCreateSubdomains2D(), PCGASMGetSubdomains() 930b1a0a8a3SJed Brown @*/ 9316a4f0f73SDmitry Karpeev PetscErrorCode PCGASMSetSubdomains(PC pc,PetscInt n,IS iis[],IS ois[]) 932b1a0a8a3SJed Brown { 933*8f3b4b4dSDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 934b1a0a8a3SJed Brown PetscErrorCode ierr; 935b1a0a8a3SJed Brown 936b1a0a8a3SJed Brown PetscFunctionBegin; 937b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 9386a4f0f73SDmitry Karpeev ierr = PetscTryMethod(pc,"PCGASMSetSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,n,iis,ois));CHKERRQ(ierr); 939*8f3b4b4dSDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 940b1a0a8a3SJed Brown PetscFunctionReturn(0); 941b1a0a8a3SJed Brown } 942b1a0a8a3SJed Brown 943b1a0a8a3SJed Brown 944b1a0a8a3SJed Brown #undef __FUNCT__ 945f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetOverlap" 946b1a0a8a3SJed Brown /*@ 947f746d493SDmitry Karpeev PCGASMSetOverlap - Sets the overlap between a pair of subdomains for the 948b1a0a8a3SJed Brown additive Schwarz preconditioner. Either all or no processors in the 949*8f3b4b4dSDmitry Karpeev pc communicator must call this routine. 950b1a0a8a3SJed Brown 951*8f3b4b4dSDmitry Karpeev Logically Collective on pc 952b1a0a8a3SJed Brown 953b1a0a8a3SJed Brown Input Parameters: 954b1a0a8a3SJed Brown + pc - the preconditioner context 955*8f3b4b4dSDmitry Karpeev - ovl - the amount of overlap between subdomains (ovl >= 0, default value = 0) 956b1a0a8a3SJed Brown 957b1a0a8a3SJed Brown Options Database Key: 95806b43e7eSDmitry Karpeev . -pc_gasm_overlap <overlap> - Sets overlap 959b1a0a8a3SJed Brown 960b1a0a8a3SJed Brown Notes: 96106b43e7eSDmitry Karpeev By default the GASM preconditioner uses 1 subdomain per processor. To use 962*8f3b4b4dSDmitry Karpeev multiple subdomain per perocessor or "straddling" subdomains that intersect 963*8f3b4b4dSDmitry Karpeev multiple processors use PCGASMSetSubdomains() (or option -pc_gasm_total_subdomains <n>). 964b1a0a8a3SJed Brown 965*8f3b4b4dSDmitry Karpeev The overlap defaults to 0, so if one desires that no additional 966b1a0a8a3SJed Brown overlap be computed beyond what may have been set with a call to 967*8f3b4b4dSDmitry Karpeev PCGASMSetSubdomains(), then ovl must be set to be 0. In particular, if one does 968*8f3b4b4dSDmitry Karpeev not explicitly set the subdomains in application code, then all overlap would be computed 969f746d493SDmitry Karpeev internally by PETSc, and using an overlap of 0 would result in an GASM 970b1a0a8a3SJed Brown variant that is equivalent to the block Jacobi preconditioner. 971b1a0a8a3SJed Brown 972b1a0a8a3SJed Brown Note that one can define initial index sets with any overlap via 97306b43e7eSDmitry Karpeev PCGASMSetSubdomains(); the routine PCGASMSetOverlap() merely allows 97406b43e7eSDmitry Karpeev PETSc to extend that overlap further, if desired. 975b1a0a8a3SJed Brown 976b1a0a8a3SJed Brown Level: intermediate 977b1a0a8a3SJed Brown 978f746d493SDmitry Karpeev .keywords: PC, GASM, set, overlap 979b1a0a8a3SJed Brown 980*8f3b4b4dSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMGetSubKSP(), 98106b43e7eSDmitry Karpeev PCGASMCreateSubdomains2D(), PCGASMGetSubdomains() 982b1a0a8a3SJed Brown @*/ 9837087cfbeSBarry Smith PetscErrorCode PCGASMSetOverlap(PC pc,PetscInt ovl) 984b1a0a8a3SJed Brown { 985b1a0a8a3SJed Brown PetscErrorCode ierr; 986*8f3b4b4dSDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 987b1a0a8a3SJed Brown 988b1a0a8a3SJed Brown PetscFunctionBegin; 989b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 990b1a0a8a3SJed Brown PetscValidLogicalCollectiveInt(pc,ovl,2); 991f746d493SDmitry Karpeev ierr = PetscTryMethod(pc,"PCGASMSetOverlap_C",(PC,PetscInt),(pc,ovl));CHKERRQ(ierr); 992*8f3b4b4dSDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 993b1a0a8a3SJed Brown PetscFunctionReturn(0); 994b1a0a8a3SJed Brown } 995b1a0a8a3SJed Brown 996b1a0a8a3SJed Brown #undef __FUNCT__ 997f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetType" 998b1a0a8a3SJed Brown /*@ 999f746d493SDmitry Karpeev PCGASMSetType - Sets the type of restriction and interpolation used 1000b1a0a8a3SJed Brown for local problems in the additive Schwarz method. 1001b1a0a8a3SJed Brown 1002b1a0a8a3SJed Brown Logically Collective on PC 1003b1a0a8a3SJed Brown 1004b1a0a8a3SJed Brown Input Parameters: 1005b1a0a8a3SJed Brown + pc - the preconditioner context 1006f746d493SDmitry Karpeev - type - variant of GASM, one of 1007b1a0a8a3SJed Brown .vb 1008f746d493SDmitry Karpeev PC_GASM_BASIC - full interpolation and restriction 1009f746d493SDmitry Karpeev PC_GASM_RESTRICT - full restriction, local processor interpolation 1010f746d493SDmitry Karpeev PC_GASM_INTERPOLATE - full interpolation, local processor restriction 1011f746d493SDmitry Karpeev PC_GASM_NONE - local processor restriction and interpolation 1012b1a0a8a3SJed Brown .ve 1013b1a0a8a3SJed Brown 1014b1a0a8a3SJed Brown Options Database Key: 1015f746d493SDmitry Karpeev . -pc_gasm_type [basic,restrict,interpolate,none] - Sets GASM type 1016b1a0a8a3SJed Brown 1017b1a0a8a3SJed Brown Level: intermediate 1018b1a0a8a3SJed Brown 1019f746d493SDmitry Karpeev .keywords: PC, GASM, set, type 1020b1a0a8a3SJed Brown 1021*8f3b4b4dSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMGetSubKSP(), 1022f746d493SDmitry Karpeev PCGASMCreateSubdomains2D() 1023b1a0a8a3SJed Brown @*/ 10247087cfbeSBarry Smith PetscErrorCode PCGASMSetType(PC pc,PCGASMType type) 1025b1a0a8a3SJed Brown { 1026b1a0a8a3SJed Brown PetscErrorCode ierr; 1027b1a0a8a3SJed Brown 1028b1a0a8a3SJed Brown PetscFunctionBegin; 1029b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1030b1a0a8a3SJed Brown PetscValidLogicalCollectiveEnum(pc,type,2); 1031f746d493SDmitry Karpeev ierr = PetscTryMethod(pc,"PCGASMSetType_C",(PC,PCGASMType),(pc,type));CHKERRQ(ierr); 1032b1a0a8a3SJed Brown PetscFunctionReturn(0); 1033b1a0a8a3SJed Brown } 1034b1a0a8a3SJed Brown 1035b1a0a8a3SJed Brown #undef __FUNCT__ 1036f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetSortIndices" 1037b1a0a8a3SJed Brown /*@ 1038f746d493SDmitry Karpeev PCGASMSetSortIndices - Determines whether subdomain indices are sorted. 1039b1a0a8a3SJed Brown 1040b1a0a8a3SJed Brown Logically Collective on PC 1041b1a0a8a3SJed Brown 1042b1a0a8a3SJed Brown Input Parameters: 1043b1a0a8a3SJed Brown + pc - the preconditioner context 1044b1a0a8a3SJed Brown - doSort - sort the subdomain indices 1045b1a0a8a3SJed Brown 1046b1a0a8a3SJed Brown Level: intermediate 1047b1a0a8a3SJed Brown 1048f746d493SDmitry Karpeev .keywords: PC, GASM, set, type 1049b1a0a8a3SJed Brown 1050*8f3b4b4dSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMGetSubKSP(), 1051f746d493SDmitry Karpeev PCGASMCreateSubdomains2D() 1052b1a0a8a3SJed Brown @*/ 10537087cfbeSBarry Smith PetscErrorCode PCGASMSetSortIndices(PC pc,PetscBool doSort) 1054b1a0a8a3SJed Brown { 1055b1a0a8a3SJed Brown PetscErrorCode ierr; 1056b1a0a8a3SJed Brown 1057b1a0a8a3SJed Brown PetscFunctionBegin; 1058b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1059b1a0a8a3SJed Brown PetscValidLogicalCollectiveBool(pc,doSort,2); 1060f746d493SDmitry Karpeev ierr = PetscTryMethod(pc,"PCGASMSetSortIndices_C",(PC,PetscBool),(pc,doSort));CHKERRQ(ierr); 1061b1a0a8a3SJed Brown PetscFunctionReturn(0); 1062b1a0a8a3SJed Brown } 1063b1a0a8a3SJed Brown 1064b1a0a8a3SJed Brown #undef __FUNCT__ 1065f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMGetSubKSP" 1066b1a0a8a3SJed Brown /*@C 1067f746d493SDmitry Karpeev PCGASMGetSubKSP - Gets the local KSP contexts for all blocks on 1068b1a0a8a3SJed Brown this processor. 1069b1a0a8a3SJed Brown 1070b1a0a8a3SJed Brown Collective on PC iff first_local is requested 1071b1a0a8a3SJed Brown 1072b1a0a8a3SJed Brown Input Parameter: 1073b1a0a8a3SJed Brown . pc - the preconditioner context 1074b1a0a8a3SJed Brown 1075b1a0a8a3SJed Brown Output Parameters: 10760298fd71SBarry Smith + n_local - the number of blocks on this processor or NULL 10770298fd71SBarry Smith . first_local - the global number of the first block on this processor or NULL, 10780298fd71SBarry Smith all processors must request or all must pass NULL 1079b1a0a8a3SJed Brown - ksp - the array of KSP contexts 1080b1a0a8a3SJed Brown 1081b1a0a8a3SJed Brown Note: 1082f746d493SDmitry Karpeev After PCGASMGetSubKSP() the array of KSPes is not to be freed 1083b1a0a8a3SJed Brown 1084b1a0a8a3SJed Brown Currently for some matrix implementations only 1 block per processor 1085b1a0a8a3SJed Brown is supported. 1086b1a0a8a3SJed Brown 1087f746d493SDmitry Karpeev You must call KSPSetUp() before calling PCGASMGetSubKSP(). 1088b1a0a8a3SJed Brown 1089b1a0a8a3SJed Brown Level: advanced 1090b1a0a8a3SJed Brown 1091f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, get, sub, KSP, context 1092b1a0a8a3SJed Brown 1093*8f3b4b4dSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMSetOverlap(), 1094f746d493SDmitry Karpeev PCGASMCreateSubdomains2D(), 1095b1a0a8a3SJed Brown @*/ 10967087cfbeSBarry Smith PetscErrorCode PCGASMGetSubKSP(PC pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp[]) 1097b1a0a8a3SJed Brown { 1098b1a0a8a3SJed Brown PetscErrorCode ierr; 1099b1a0a8a3SJed Brown 1100b1a0a8a3SJed Brown PetscFunctionBegin; 1101b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1102f746d493SDmitry Karpeev ierr = PetscUseMethod(pc,"PCGASMGetSubKSP_C",(PC,PetscInt*,PetscInt*,KSP **),(pc,n_local,first_local,ksp));CHKERRQ(ierr); 1103b1a0a8a3SJed Brown PetscFunctionReturn(0); 1104b1a0a8a3SJed Brown } 1105b1a0a8a3SJed Brown 1106b1a0a8a3SJed Brown /* -------------------------------------------------------------------------------------*/ 1107b1a0a8a3SJed Brown /*MC 1108f746d493SDmitry Karpeev PCGASM - Use the (restricted) additive Schwarz method, each block is (approximately) solved with 1109b1a0a8a3SJed Brown its own KSP object. 1110b1a0a8a3SJed Brown 1111b1a0a8a3SJed Brown Options Database Keys: 1112*8f3b4b4dSDmitry Karpeev + -pc_gasm_total_subdomains <n> - Sets total number of local subdomains to be distributed among processors 111306b43e7eSDmitry Karpeev . -pc_gasm_view_subdomains - activates the printing of subdomain indices in PCView(), -ksp_view or -snes_view 111406b43e7eSDmitry Karpeev . -pc_gasm_print_subdomains - activates the printing of subdomain indices in PCSetUp() 111506b43e7eSDmitry Karpeev . -pc_gasm_overlap <ovl> - Sets overlap by which to (automatically) extend local subdomains 1116f746d493SDmitry Karpeev - -pc_gasm_type [basic,restrict,interpolate,none] - Sets GASM type 1117b1a0a8a3SJed Brown 1118b1a0a8a3SJed Brown IMPORTANT: If you run with, for example, 3 blocks on 1 processor or 3 blocks on 3 processors you 1119f746d493SDmitry Karpeev will get a different convergence rate due to the default option of -pc_gasm_type restrict. Use 1120f746d493SDmitry Karpeev -pc_gasm_type basic to use the standard GASM. 1121b1a0a8a3SJed Brown 1122b1a0a8a3SJed Brown Notes: Each processor can have one or more blocks, but a block cannot be shared by more 1123b1a0a8a3SJed Brown than one processor. Defaults to one block per processor. 1124b1a0a8a3SJed Brown 1125b1a0a8a3SJed Brown To set options on the solvers for each block append -sub_ to all the KSP, and PC 1126b1a0a8a3SJed Brown options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly 1127b1a0a8a3SJed Brown 1128f746d493SDmitry Karpeev To set the options on the solvers separate for each block call PCGASMGetSubKSP() 1129b1a0a8a3SJed Brown and set the options directly on the resulting KSP object (you can access its PC 1130b1a0a8a3SJed Brown with KSPGetPC()) 1131b1a0a8a3SJed Brown 1132b1a0a8a3SJed Brown 1133b1a0a8a3SJed Brown Level: beginner 1134b1a0a8a3SJed Brown 1135b1a0a8a3SJed Brown Concepts: additive Schwarz method 1136b1a0a8a3SJed Brown 1137b1a0a8a3SJed Brown References: 1138b1a0a8a3SJed Brown An additive variant of the Schwarz alternating method for the case of many subregions 1139b1a0a8a3SJed Brown M Dryja, OB Widlund - Courant Institute, New York University Technical report 1140b1a0a8a3SJed Brown 1141b1a0a8a3SJed Brown Domain Decompositions: Parallel Multilevel Methods for Elliptic Partial Differential Equations, 1142b1a0a8a3SJed Brown Barry Smith, Petter Bjorstad, and William Gropp, Cambridge University Press, ISBN 0-521-49589-X. 1143b1a0a8a3SJed Brown 1144b1a0a8a3SJed Brown .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 114549517cdeSBarry Smith PCBJACOBI, PCGASMGetSubKSP(), PCGASMSetSubdomains(), 1146*8f3b4b4dSDmitry Karpeev PCSetModifySubmatrices(), PCGASMSetOverlap(), PCGASMSetType() 1147b1a0a8a3SJed Brown 1148b1a0a8a3SJed Brown M*/ 1149b1a0a8a3SJed Brown 1150b1a0a8a3SJed Brown #undef __FUNCT__ 1151f746d493SDmitry Karpeev #define __FUNCT__ "PCCreate_GASM" 11528cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_GASM(PC pc) 1153b1a0a8a3SJed Brown { 1154b1a0a8a3SJed Brown PetscErrorCode ierr; 1155f746d493SDmitry Karpeev PC_GASM *osm; 1156b1a0a8a3SJed Brown 1157b1a0a8a3SJed Brown PetscFunctionBegin; 1158b00a9115SJed Brown ierr = PetscNewLog(pc,&osm);CHKERRQ(ierr); 11592fa5cd67SKarl Rupp 1160*8f3b4b4dSDmitry Karpeev osm->N = PETSC_DETERMINE; 116106b43e7eSDmitry Karpeev osm->n = PETSC_DECIDE; 1162*8f3b4b4dSDmitry Karpeev osm->nmax = PETSC_DETERMINE; 1163*8f3b4b4dSDmitry Karpeev osm->overlap = 0; 1164b1a0a8a3SJed Brown osm->ksp = 0; 11656a4f0f73SDmitry Karpeev osm->gorestriction = 0; 11666a4f0f73SDmitry Karpeev osm->girestriction = 0; 1167ab3e923fSDmitry Karpeev osm->gx = 0; 1168ab3e923fSDmitry Karpeev osm->gy = 0; 1169b1a0a8a3SJed Brown osm->x = 0; 1170b1a0a8a3SJed Brown osm->y = 0; 11716a4f0f73SDmitry Karpeev osm->ois = 0; 11726a4f0f73SDmitry Karpeev osm->iis = 0; 1173b1a0a8a3SJed Brown osm->pmat = 0; 1174f746d493SDmitry Karpeev osm->type = PC_GASM_RESTRICT; 11756a4f0f73SDmitry Karpeev osm->same_subdomain_solvers = PETSC_TRUE; 1176b1a0a8a3SJed Brown osm->sort_indices = PETSC_TRUE; 1177d709ea83SDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 1178b1a0a8a3SJed Brown 1179b1a0a8a3SJed Brown pc->data = (void*)osm; 1180f746d493SDmitry Karpeev pc->ops->apply = PCApply_GASM; 1181f746d493SDmitry Karpeev pc->ops->applytranspose = PCApplyTranspose_GASM; 1182f746d493SDmitry Karpeev pc->ops->setup = PCSetUp_GASM; 1183a06653b4SBarry Smith pc->ops->reset = PCReset_GASM; 1184f746d493SDmitry Karpeev pc->ops->destroy = PCDestroy_GASM; 1185f746d493SDmitry Karpeev pc->ops->setfromoptions = PCSetFromOptions_GASM; 1186f746d493SDmitry Karpeev pc->ops->setuponblocks = PCSetUpOnBlocks_GASM; 1187f746d493SDmitry Karpeev pc->ops->view = PCView_GASM; 1188b1a0a8a3SJed Brown pc->ops->applyrichardson = 0; 1189b1a0a8a3SJed Brown 1190bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetSubdomains_C",PCGASMSetSubdomains_GASM);CHKERRQ(ierr); 1191bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetOverlap_C",PCGASMSetOverlap_GASM);CHKERRQ(ierr); 1192bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetType_C",PCGASMSetType_GASM);CHKERRQ(ierr); 1193bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetSortIndices_C",PCGASMSetSortIndices_GASM);CHKERRQ(ierr); 1194bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGASMGetSubKSP_C",PCGASMGetSubKSP_GASM);CHKERRQ(ierr); 1195b1a0a8a3SJed Brown PetscFunctionReturn(0); 1196b1a0a8a3SJed Brown } 1197b1a0a8a3SJed Brown 1198b1a0a8a3SJed Brown 1199b1a0a8a3SJed Brown #undef __FUNCT__ 120006b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMCreateLocalSubdomains" 1201*8f3b4b4dSDmitry Karpeev PetscErrorCode PCGASMCreateLocalSubdomains(Mat A, PetscInt nloc, IS *iis[]) 1202b1a0a8a3SJed Brown { 1203b1a0a8a3SJed Brown MatPartitioning mpart; 1204b1a0a8a3SJed Brown const char *prefix; 120511bd1e4dSLisandro Dalcin PetscErrorCode (*f)(Mat,MatReuse,Mat*); 1206b1a0a8a3SJed Brown PetscMPIInt size; 1207b1a0a8a3SJed Brown PetscInt i,j,rstart,rend,bs; 120811bd1e4dSLisandro Dalcin PetscBool isbaij = PETSC_FALSE,foundpart = PETSC_FALSE; 12090298fd71SBarry Smith Mat Ad = NULL, adj; 1210b1a0a8a3SJed Brown IS ispart,isnumb,*is; 1211b1a0a8a3SJed Brown PetscErrorCode ierr; 1212b1a0a8a3SJed Brown 1213b1a0a8a3SJed Brown PetscFunctionBegin; 1214*8f3b4b4dSDmitry Karpeev if (nloc < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of local subdomains must > 0, got nloc = %D",nloc); 1215b1a0a8a3SJed Brown 1216b1a0a8a3SJed Brown /* Get prefix, row distribution, and block size */ 1217b1a0a8a3SJed Brown ierr = MatGetOptionsPrefix(A,&prefix);CHKERRQ(ierr); 1218b1a0a8a3SJed Brown ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 1219b1a0a8a3SJed Brown ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr); 1220b1a0a8a3SJed 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); 1221b1a0a8a3SJed Brown 1222b1a0a8a3SJed Brown /* Get diagonal block from matrix if possible */ 1223ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr); 12240005d66cSJed Brown ierr = PetscObjectQueryFunction((PetscObject)A,"MatGetDiagonalBlock_C",&f);CHKERRQ(ierr); 1225b1a0a8a3SJed Brown if (f) { 122611bd1e4dSLisandro Dalcin ierr = MatGetDiagonalBlock(A,&Ad);CHKERRQ(ierr); 1227b1a0a8a3SJed Brown } else if (size == 1) { 122811bd1e4dSLisandro Dalcin Ad = A; 1229b1a0a8a3SJed Brown } 1230b1a0a8a3SJed Brown if (Ad) { 1231251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)Ad,MATSEQBAIJ,&isbaij);CHKERRQ(ierr); 1232251f4c67SDmitry Karpeev if (!isbaij) {ierr = PetscObjectTypeCompare((PetscObject)Ad,MATSEQSBAIJ,&isbaij);CHKERRQ(ierr);} 1233b1a0a8a3SJed Brown } 1234*8f3b4b4dSDmitry Karpeev if (Ad && nloc > 1) { 1235b1a0a8a3SJed Brown PetscBool match,done; 1236b1a0a8a3SJed Brown /* Try to setup a good matrix partitioning if available */ 1237b1a0a8a3SJed Brown ierr = MatPartitioningCreate(PETSC_COMM_SELF,&mpart);CHKERRQ(ierr); 1238b1a0a8a3SJed Brown ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr); 1239b1a0a8a3SJed Brown ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr); 1240251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGCURRENT,&match);CHKERRQ(ierr); 1241b1a0a8a3SJed Brown if (!match) { 1242251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGSQUARE,&match);CHKERRQ(ierr); 1243b1a0a8a3SJed Brown } 1244b1a0a8a3SJed Brown if (!match) { /* assume a "good" partitioner is available */ 12451a83f524SJed Brown PetscInt na; 12461a83f524SJed Brown const PetscInt *ia,*ja; 1247b1a0a8a3SJed Brown ierr = MatGetRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr); 1248b1a0a8a3SJed Brown if (done) { 1249b1a0a8a3SJed Brown /* Build adjacency matrix by hand. Unfortunately a call to 1250b1a0a8a3SJed Brown MatConvert(Ad,MATMPIADJ,MAT_INITIAL_MATRIX,&adj) will 1251b1a0a8a3SJed Brown remove the block-aij structure and we cannot expect 1252b1a0a8a3SJed Brown MatPartitioning to split vertices as we need */ 12531a83f524SJed Brown PetscInt i,j,len,nnz,cnt,*iia=0,*jja=0; 12541a83f524SJed Brown const PetscInt *row; 1255b1a0a8a3SJed Brown nnz = 0; 1256b1a0a8a3SJed Brown for (i=0; i<na; i++) { /* count number of nonzeros */ 1257b1a0a8a3SJed Brown len = ia[i+1] - ia[i]; 1258b1a0a8a3SJed Brown row = ja + ia[i]; 1259b1a0a8a3SJed Brown for (j=0; j<len; j++) { 1260b1a0a8a3SJed Brown if (row[j] == i) { /* don't count diagonal */ 1261b1a0a8a3SJed Brown len--; break; 1262b1a0a8a3SJed Brown } 1263b1a0a8a3SJed Brown } 1264b1a0a8a3SJed Brown nnz += len; 1265b1a0a8a3SJed Brown } 1266854ce69bSBarry Smith ierr = PetscMalloc1(na+1,&iia);CHKERRQ(ierr); 1267854ce69bSBarry Smith ierr = PetscMalloc1(nnz,&jja);CHKERRQ(ierr); 1268b1a0a8a3SJed Brown nnz = 0; 1269b1a0a8a3SJed Brown iia[0] = 0; 1270b1a0a8a3SJed Brown for (i=0; i<na; i++) { /* fill adjacency */ 1271b1a0a8a3SJed Brown cnt = 0; 1272b1a0a8a3SJed Brown len = ia[i+1] - ia[i]; 1273b1a0a8a3SJed Brown row = ja + ia[i]; 1274b1a0a8a3SJed Brown for (j=0; j<len; j++) { 12752fa5cd67SKarl Rupp if (row[j] != i) jja[nnz+cnt++] = row[j]; /* if not diagonal */ 1276b1a0a8a3SJed Brown } 1277b1a0a8a3SJed Brown nnz += cnt; 1278b1a0a8a3SJed Brown iia[i+1] = nnz; 1279b1a0a8a3SJed Brown } 1280b1a0a8a3SJed Brown /* Partitioning of the adjacency matrix */ 12810298fd71SBarry Smith ierr = MatCreateMPIAdj(PETSC_COMM_SELF,na,na,iia,jja,NULL,&adj);CHKERRQ(ierr); 1282b1a0a8a3SJed Brown ierr = MatPartitioningSetAdjacency(mpart,adj);CHKERRQ(ierr); 1283*8f3b4b4dSDmitry Karpeev ierr = MatPartitioningSetNParts(mpart,nloc);CHKERRQ(ierr); 1284b1a0a8a3SJed Brown ierr = MatPartitioningApply(mpart,&ispart);CHKERRQ(ierr); 1285b1a0a8a3SJed Brown ierr = ISPartitioningToNumbering(ispart,&isnumb);CHKERRQ(ierr); 12866bf464f9SBarry Smith ierr = MatDestroy(&adj);CHKERRQ(ierr); 1287b1a0a8a3SJed Brown foundpart = PETSC_TRUE; 1288b1a0a8a3SJed Brown } 1289b1a0a8a3SJed Brown ierr = MatRestoreRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr); 1290b1a0a8a3SJed Brown } 12916bf464f9SBarry Smith ierr = MatPartitioningDestroy(&mpart);CHKERRQ(ierr); 1292b1a0a8a3SJed Brown } 1293*8f3b4b4dSDmitry Karpeev ierr = PetscMalloc1(nloc,&is);CHKERRQ(ierr); 1294b1a0a8a3SJed Brown if (!foundpart) { 1295b1a0a8a3SJed Brown 1296b1a0a8a3SJed Brown /* Partitioning by contiguous chunks of rows */ 1297b1a0a8a3SJed Brown 1298b1a0a8a3SJed Brown PetscInt mbs = (rend-rstart)/bs; 1299b1a0a8a3SJed Brown PetscInt start = rstart; 1300*8f3b4b4dSDmitry Karpeev for (i=0; i<nloc; i++) { 1301*8f3b4b4dSDmitry Karpeev PetscInt count = (mbs/nloc + ((mbs % nloc) > i)) * bs; 1302b1a0a8a3SJed Brown ierr = ISCreateStride(PETSC_COMM_SELF,count,start,1,&is[i]);CHKERRQ(ierr); 1303b1a0a8a3SJed Brown start += count; 1304b1a0a8a3SJed Brown } 1305b1a0a8a3SJed Brown 1306b1a0a8a3SJed Brown } else { 1307b1a0a8a3SJed Brown 1308b1a0a8a3SJed Brown /* Partitioning by adjacency of diagonal block */ 1309b1a0a8a3SJed Brown 1310b1a0a8a3SJed Brown const PetscInt *numbering; 1311b1a0a8a3SJed Brown PetscInt *count,nidx,*indices,*newidx,start=0; 1312b1a0a8a3SJed Brown /* Get node count in each partition */ 1313*8f3b4b4dSDmitry Karpeev ierr = PetscMalloc1(nloc,&count);CHKERRQ(ierr); 1314*8f3b4b4dSDmitry Karpeev ierr = ISPartitioningCount(ispart,nloc,count);CHKERRQ(ierr); 1315b1a0a8a3SJed Brown if (isbaij && bs > 1) { /* adjust for the block-aij case */ 1316*8f3b4b4dSDmitry Karpeev for (i=0; i<nloc; i++) count[i] *= bs; 1317b1a0a8a3SJed Brown } 1318b1a0a8a3SJed Brown /* Build indices from node numbering */ 1319b1a0a8a3SJed Brown ierr = ISGetLocalSize(isnumb,&nidx);CHKERRQ(ierr); 1320785e854fSJed Brown ierr = PetscMalloc1(nidx,&indices);CHKERRQ(ierr); 1321b1a0a8a3SJed Brown for (i=0; i<nidx; i++) indices[i] = i; /* needs to be initialized */ 1322b1a0a8a3SJed Brown ierr = ISGetIndices(isnumb,&numbering);CHKERRQ(ierr); 1323b1a0a8a3SJed Brown ierr = PetscSortIntWithPermutation(nidx,numbering,indices);CHKERRQ(ierr); 1324b1a0a8a3SJed Brown ierr = ISRestoreIndices(isnumb,&numbering);CHKERRQ(ierr); 1325b1a0a8a3SJed Brown if (isbaij && bs > 1) { /* adjust for the block-aij case */ 1326785e854fSJed Brown ierr = PetscMalloc1(nidx*bs,&newidx);CHKERRQ(ierr); 13272fa5cd67SKarl Rupp for (i=0; i<nidx; i++) { 13282fa5cd67SKarl Rupp for (j=0; j<bs; j++) newidx[i*bs+j] = indices[i]*bs + j; 13292fa5cd67SKarl Rupp } 1330b1a0a8a3SJed Brown ierr = PetscFree(indices);CHKERRQ(ierr); 1331b1a0a8a3SJed Brown nidx *= bs; 1332b1a0a8a3SJed Brown indices = newidx; 1333b1a0a8a3SJed Brown } 1334b1a0a8a3SJed Brown /* Shift to get global indices */ 1335b1a0a8a3SJed Brown for (i=0; i<nidx; i++) indices[i] += rstart; 1336b1a0a8a3SJed Brown 1337b1a0a8a3SJed Brown /* Build the index sets for each block */ 1338*8f3b4b4dSDmitry Karpeev for (i=0; i<nloc; i++) { 1339b1a0a8a3SJed Brown ierr = ISCreateGeneral(PETSC_COMM_SELF,count[i],&indices[start],PETSC_COPY_VALUES,&is[i]);CHKERRQ(ierr); 1340b1a0a8a3SJed Brown ierr = ISSort(is[i]);CHKERRQ(ierr); 1341b1a0a8a3SJed Brown start += count[i]; 1342b1a0a8a3SJed Brown } 1343b1a0a8a3SJed Brown 1344302440fdSBarry Smith ierr = PetscFree(count);CHKERRQ(ierr); 1345302440fdSBarry Smith ierr = PetscFree(indices);CHKERRQ(ierr); 1346fcfd50ebSBarry Smith ierr = ISDestroy(&isnumb);CHKERRQ(ierr); 1347fcfd50ebSBarry Smith ierr = ISDestroy(&ispart);CHKERRQ(ierr); 1348b1a0a8a3SJed Brown } 13496a4f0f73SDmitry Karpeev *iis = is; 1350*8f3b4b4dSDmitry Karpeev PetscFunctionReturn(0); 1351*8f3b4b4dSDmitry Karpeev } 1352*8f3b4b4dSDmitry Karpeev 1353*8f3b4b4dSDmitry Karpeev #undef __FUNCT__ 1354*8f3b4b4dSDmitry Karpeev #define __FUNCT__ "PCGASMCreateStraddlingSubdomains" 1355*8f3b4b4dSDmitry Karpeev PETSC_INTERN PetscErrorCode PCGASMCreateStraddlingSubdomains(Mat A,PetscMPIInt N,PetscMPIInt *n,IS *iis[]) 1356*8f3b4b4dSDmitry Karpeev { 1357*8f3b4b4dSDmitry Karpeev PetscErrorCode ierr; 1358*8f3b4b4dSDmitry Karpeev 1359*8f3b4b4dSDmitry Karpeev PetscFunctionBegin; 1360*8f3b4b4dSDmitry Karpeev ierr = MatSubdomainsCreateCoalesce(A,N,n,iis);CHKERRQ(ierr); 1361*8f3b4b4dSDmitry Karpeev PetscFunctionReturn(0); 1362*8f3b4b4dSDmitry Karpeev } 1363*8f3b4b4dSDmitry Karpeev 1364*8f3b4b4dSDmitry Karpeev 1365*8f3b4b4dSDmitry Karpeev 1366*8f3b4b4dSDmitry Karpeev #undef __FUNCT__ 1367*8f3b4b4dSDmitry Karpeev #define __FUNCT__ "PCGASMCreateSubdomains" 1368*8f3b4b4dSDmitry Karpeev /*@C 1369*8f3b4b4dSDmitry Karpeev PCGASMCreateSubdomains - Creates n index sets defining n nonoverlapping subdomains for the additive 1370*8f3b4b4dSDmitry Karpeev Schwarz preconditioner for a any problem based on its matrix. 1371*8f3b4b4dSDmitry Karpeev 1372*8f3b4b4dSDmitry Karpeev Collective 1373*8f3b4b4dSDmitry Karpeev 1374*8f3b4b4dSDmitry Karpeev Input Parameters: 1375*8f3b4b4dSDmitry Karpeev + A - The global matrix operator 1376*8f3b4b4dSDmitry Karpeev - N - the number of global subdomains requested 1377*8f3b4b4dSDmitry Karpeev 1378*8f3b4b4dSDmitry Karpeev Output Parameters: 1379*8f3b4b4dSDmitry Karpeev + n - the number of subdomains created on this processor 1380*8f3b4b4dSDmitry Karpeev - iis - the array of index sets defining the local inner subdomains (on which the correction is applied) 1381*8f3b4b4dSDmitry Karpeev 1382*8f3b4b4dSDmitry Karpeev Level: advanced 1383*8f3b4b4dSDmitry Karpeev 1384*8f3b4b4dSDmitry Karpeev Note: When N >= A's communicator size, each subdomain is local -- contained within a single processor. 1385*8f3b4b4dSDmitry Karpeev When N < size, the subdomains are 'straddling' (processor boundaries) and are no longer local. 1386*8f3b4b4dSDmitry Karpeev The resulting subdomains can be use in PCGASMSetSubdomains(pc,n,iss,NULL). The overlapping 1387*8f3b4b4dSDmitry Karpeev outer subdomains will be automatically generated from these according to the requested amount of 1388*8f3b4b4dSDmitry Karpeev overlap; this is currently supported only with local subdomains. 1389*8f3b4b4dSDmitry Karpeev 1390*8f3b4b4dSDmitry Karpeev 1391*8f3b4b4dSDmitry Karpeev .keywords: PC, GASM, additive Schwarz, create, subdomains, unstructured grid 1392*8f3b4b4dSDmitry Karpeev 1393*8f3b4b4dSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMDestroySubdomains() 1394*8f3b4b4dSDmitry Karpeev @*/ 1395*8f3b4b4dSDmitry Karpeev PetscErrorCode PCGASMCreateSubdomains(Mat A, PetscMPIInt N, PetscMPIInt *n, IS *iis[]) 1396*8f3b4b4dSDmitry Karpeev { 1397*8f3b4b4dSDmitry Karpeev PetscMPIInt size; 1398*8f3b4b4dSDmitry Karpeev PetscErrorCode ierr; 1399*8f3b4b4dSDmitry Karpeev 1400*8f3b4b4dSDmitry Karpeev PetscFunctionBegin; 1401*8f3b4b4dSDmitry Karpeev PetscValidHeaderSpecific(A,MAT_CLASSID,1); 1402*8f3b4b4dSDmitry Karpeev PetscValidPointer(iis,4); 1403*8f3b4b4dSDmitry Karpeev 1404*8f3b4b4dSDmitry Karpeev if (N < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Number of subdomains must be > 0, N = %D",N); 1405*8f3b4b4dSDmitry Karpeev ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr); 1406*8f3b4b4dSDmitry Karpeev if (N >= size) { 1407*8f3b4b4dSDmitry Karpeev *n = N/size + (N%size); 1408*8f3b4b4dSDmitry Karpeev ierr = PCGASMCreateLocalSubdomains(A,*n,iis);CHKERRQ(ierr); 14096a4f0f73SDmitry Karpeev } else { 1410*8f3b4b4dSDmitry Karpeev ierr = PCGASMCreateStraddlingSubdomains(A,N,n,iis);CHKERRQ(ierr); 14116a4f0f73SDmitry Karpeev } 1412b1a0a8a3SJed Brown PetscFunctionReturn(0); 1413b1a0a8a3SJed Brown } 1414b1a0a8a3SJed Brown 1415b1a0a8a3SJed Brown #undef __FUNCT__ 1416f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMDestroySubdomains" 1417b1a0a8a3SJed Brown /*@C 1418f746d493SDmitry Karpeev PCGASMDestroySubdomains - Destroys the index sets created with 1419*8f3b4b4dSDmitry Karpeev PCGASMCreateSubdomains() or PCGASMCreateSubdomains2D. Should be 142006b43e7eSDmitry Karpeev called after setting subdomains with PCGASMSetSubdomains(). 1421b1a0a8a3SJed Brown 1422b1a0a8a3SJed Brown Collective 1423b1a0a8a3SJed Brown 1424b1a0a8a3SJed Brown Input Parameters: 1425b1a0a8a3SJed Brown + n - the number of index sets 14266a4f0f73SDmitry Karpeev . iis - the array of inner subdomains, 14270298fd71SBarry Smith - ois - the array of outer subdomains, can be NULL 1428b1a0a8a3SJed Brown 14296a4f0f73SDmitry Karpeev Level: intermediate 14306a4f0f73SDmitry Karpeev 14316a4f0f73SDmitry Karpeev Notes: this is merely a convenience subroutine that walks each list, 14326a4f0f73SDmitry Karpeev destroys each IS on the list, and then frees the list. 1433b1a0a8a3SJed Brown 1434f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, create, subdomains, unstructured grid 1435b1a0a8a3SJed Brown 1436*8f3b4b4dSDmitry Karpeev .seealso: PCGASMCreateSubdomains(), PCGASMSetSubdomains() 1437b1a0a8a3SJed Brown @*/ 14386a4f0f73SDmitry Karpeev PetscErrorCode PCGASMDestroySubdomains(PetscInt n, IS iis[], IS ois[]) 1439b1a0a8a3SJed Brown { 1440b1a0a8a3SJed Brown PetscInt i; 1441b1a0a8a3SJed Brown PetscErrorCode ierr; 14425fd66863SKarl Rupp 1443b1a0a8a3SJed Brown PetscFunctionBegin; 1444a35b7b57SDmitry Karpeev if (n <= 0) PetscFunctionReturn(0); 1445a35b7b57SDmitry Karpeev if (iis) { 14466a4f0f73SDmitry Karpeev PetscValidPointer(iis,2); 1447a35b7b57SDmitry Karpeev for (i=0; i<n; i++) { 1448a35b7b57SDmitry Karpeev ierr = ISDestroy(&iis[i]);CHKERRQ(ierr); 1449a35b7b57SDmitry Karpeev } 14506a4f0f73SDmitry Karpeev ierr = PetscFree(iis);CHKERRQ(ierr); 1451a35b7b57SDmitry Karpeev } 14526a4f0f73SDmitry Karpeev if (ois) { 1453a35b7b57SDmitry Karpeev for (i=0; i<n; i++) { 1454a35b7b57SDmitry Karpeev ierr = ISDestroy(&ois[i]);CHKERRQ(ierr); 1455a35b7b57SDmitry Karpeev } 14566a4f0f73SDmitry Karpeev ierr = PetscFree(ois);CHKERRQ(ierr); 1457b1a0a8a3SJed Brown } 1458b1a0a8a3SJed Brown PetscFunctionReturn(0); 1459b1a0a8a3SJed Brown } 1460b1a0a8a3SJed Brown 1461af538404SDmitry Karpeev 1462af538404SDmitry Karpeev #define PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,xleft_loc,ylow_loc,xright_loc,yhigh_loc,n) \ 1463af538404SDmitry Karpeev { \ 1464af538404SDmitry Karpeev PetscInt first_row = first/M, last_row = last/M+1; \ 1465af538404SDmitry Karpeev /* \ 1466af538404SDmitry Karpeev Compute ylow_loc and yhigh_loc so that (ylow_loc,xleft) and (yhigh_loc,xright) are the corners \ 1467af538404SDmitry Karpeev of the bounding box of the intersection of the subdomain with the local ownership range (local \ 1468af538404SDmitry Karpeev subdomain). \ 1469af538404SDmitry Karpeev Also compute xleft_loc and xright_loc as the lower and upper bounds on the first and last rows \ 1470af538404SDmitry Karpeev of the intersection. \ 1471af538404SDmitry Karpeev */ \ 1472af538404SDmitry Karpeev /* ylow_loc is the grid row containing the first element of the local sumbdomain */ \ 1473eec7e3faSDmitry Karpeev *ylow_loc = PetscMax(first_row,ylow); \ 1474af538404SDmitry Karpeev /* xleft_loc is the offset of first element of the local subdomain within its grid row (might actually be outside the local subdomain) */ \ 1475eec7e3faSDmitry Karpeev *xleft_loc = *ylow_loc==first_row ? PetscMax(first%M,xleft) : xleft; \ 1476af538404SDmitry Karpeev /* yhigh_loc is the grid row above the last local subdomain element */ \ 1477eec7e3faSDmitry Karpeev *yhigh_loc = PetscMin(last_row,yhigh); \ 1478af538404SDmitry Karpeev /* xright is the offset of the end of the local subdomain within its grid row (might actually be outside the local subdomain) */ \ 1479eec7e3faSDmitry Karpeev *xright_loc = *yhigh_loc==last_row ? PetscMin(xright,last%M) : xright; \ 1480af538404SDmitry Karpeev /* Now compute the size of the local subdomain n. */ \ 1481c3518bceSDmitry Karpeev *n = 0; \ 1482eec7e3faSDmitry Karpeev if (*ylow_loc < *yhigh_loc) { \ 1483af538404SDmitry Karpeev PetscInt width = xright-xleft; \ 1484eec7e3faSDmitry Karpeev *n += width*(*yhigh_loc-*ylow_loc-1); \ 1485eec7e3faSDmitry Karpeev *n += PetscMin(PetscMax(*xright_loc-xleft,0),width); \ 1486eec7e3faSDmitry Karpeev *n -= PetscMin(PetscMax(*xleft_loc-xleft,0), width); \ 1487af538404SDmitry Karpeev } \ 1488af538404SDmitry Karpeev } 1489af538404SDmitry Karpeev 1490c3518bceSDmitry Karpeev 1491eec7e3faSDmitry Karpeev 1492eec7e3faSDmitry Karpeev #undef __FUNCT__ 1493f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMCreateSubdomains2D" 1494b1a0a8a3SJed Brown /*@ 1495f746d493SDmitry Karpeev PCGASMCreateSubdomains2D - Creates the index sets for the overlapping Schwarz 1496b1a0a8a3SJed Brown preconditioner for a two-dimensional problem on a regular grid. 1497b1a0a8a3SJed Brown 1498af538404SDmitry Karpeev Collective 1499b1a0a8a3SJed Brown 1500b1a0a8a3SJed Brown Input Parameters: 150106b43e7eSDmitry Karpeev + M, N - the global number of grid points in the x and y directions 1502af538404SDmitry Karpeev . Mdomains, Ndomains - the global number of subdomains in the x and y directions 1503b1a0a8a3SJed Brown . dof - degrees of freedom per node 1504b1a0a8a3SJed Brown - overlap - overlap in mesh lines 1505b1a0a8a3SJed Brown 1506b1a0a8a3SJed Brown Output Parameters: 1507af538404SDmitry Karpeev + Nsub - the number of local subdomains created 15086a4f0f73SDmitry Karpeev . iis - array of index sets defining inner (nonoverlapping) subdomains 15096a4f0f73SDmitry Karpeev - ois - array of index sets defining outer (overlapping, if overlap > 0) subdomains 1510b1a0a8a3SJed Brown 1511b1a0a8a3SJed Brown 1512b1a0a8a3SJed Brown Level: advanced 1513b1a0a8a3SJed Brown 1514f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, create, subdomains, 2D, regular grid 1515b1a0a8a3SJed Brown 1516*8f3b4b4dSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMGetSubKSP(), PCGASMSetOverlap() 1517b1a0a8a3SJed Brown @*/ 15186a4f0f73SDmitry Karpeev PetscErrorCode PCGASMCreateSubdomains2D(PC pc, PetscInt M,PetscInt N,PetscInt Mdomains,PetscInt Ndomains,PetscInt dof,PetscInt overlap, PetscInt *nsub,IS **iis,IS **ois) 1519b1a0a8a3SJed Brown { 1520b1a0a8a3SJed Brown PetscErrorCode ierr; 15212bbb417fSDmitry Karpeev PetscMPIInt size, rank; 15222bbb417fSDmitry Karpeev PetscInt i, j; 15232bbb417fSDmitry Karpeev PetscInt maxheight, maxwidth; 15242bbb417fSDmitry Karpeev PetscInt xstart, xleft, xright, xleft_loc, xright_loc; 15252bbb417fSDmitry Karpeev PetscInt ystart, ylow, yhigh, ylow_loc, yhigh_loc; 1526eec7e3faSDmitry Karpeev PetscInt x[2][2], y[2][2], n[2]; 15272bbb417fSDmitry Karpeev PetscInt first, last; 15282bbb417fSDmitry Karpeev PetscInt nidx, *idx; 15292bbb417fSDmitry Karpeev PetscInt ii,jj,s,q,d; 1530af538404SDmitry Karpeev PetscInt k,kk; 15312bbb417fSDmitry Karpeev PetscMPIInt color; 15322bbb417fSDmitry Karpeev MPI_Comm comm, subcomm; 15336a4f0f73SDmitry Karpeev IS **xis = 0, **is = ois, **is_local = iis; 1534d34fcf5fSBarry Smith 1535b1a0a8a3SJed Brown PetscFunctionBegin; 15362bbb417fSDmitry Karpeev ierr = PetscObjectGetComm((PetscObject)pc, &comm);CHKERRQ(ierr); 15372bbb417fSDmitry Karpeev ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 15382bbb417fSDmitry Karpeev ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 15392bbb417fSDmitry Karpeev ierr = MatGetOwnershipRange(pc->pmat, &first, &last);CHKERRQ(ierr); 1540d34fcf5fSBarry 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) " 15412bbb417fSDmitry Karpeev "does not respect the number of degrees of freedom per grid point %D", first, last, dof); 1542d34fcf5fSBarry Smith 1543af538404SDmitry Karpeev /* Determine the number of domains with nonzero intersections with the local ownership range. */ 15442bbb417fSDmitry Karpeev s = 0; 1545af538404SDmitry Karpeev ystart = 0; 1546af538404SDmitry Karpeev for (j=0; j<Ndomains; ++j) { 1547af538404SDmitry Karpeev maxheight = N/Ndomains + ((N % Ndomains) > j); /* Maximal height of subdomain */ 1548af538404SDmitry 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); 1549eec7e3faSDmitry Karpeev /* Vertical domain limits with an overlap. */ 1550eec7e3faSDmitry Karpeev ylow = PetscMax(ystart - overlap,0); 1551eec7e3faSDmitry Karpeev yhigh = PetscMin(ystart + maxheight + overlap,N); 1552b1a0a8a3SJed Brown xstart = 0; 1553af538404SDmitry Karpeev for (i=0; i<Mdomains; ++i) { 1554af538404SDmitry Karpeev maxwidth = M/Mdomains + ((M % Mdomains) > i); /* Maximal width of subdomain */ 1555af538404SDmitry 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); 1556eec7e3faSDmitry Karpeev /* Horizontal domain limits with an overlap. */ 1557eec7e3faSDmitry Karpeev xleft = PetscMax(xstart - overlap,0); 1558eec7e3faSDmitry Karpeev xright = PetscMin(xstart + maxwidth + overlap,M); 1559af538404SDmitry Karpeev /* 1560af538404SDmitry Karpeev Determine whether this subdomain intersects this processor's ownership range of pc->pmat. 1561af538404SDmitry Karpeev */ 1562c3518bceSDmitry Karpeev PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,(&xleft_loc),(&ylow_loc),(&xright_loc),(&yhigh_loc),(&nidx)); 15632fa5cd67SKarl Rupp if (nidx) ++s; 1564af538404SDmitry Karpeev xstart += maxwidth; 1565af538404SDmitry Karpeev } /* for (i = 0; i < Mdomains; ++i) */ 1566af538404SDmitry Karpeev ystart += maxheight; 1567af538404SDmitry Karpeev } /* for (j = 0; j < Ndomains; ++j) */ 15682fa5cd67SKarl Rupp 1569af538404SDmitry Karpeev /* Now we can allocate the necessary number of ISs. */ 1570af538404SDmitry Karpeev *nsub = s; 1571854ce69bSBarry Smith ierr = PetscMalloc1(*nsub,is);CHKERRQ(ierr); 1572854ce69bSBarry Smith ierr = PetscMalloc1(*nsub,is_local);CHKERRQ(ierr); 1573af538404SDmitry Karpeev s = 0; 1574af538404SDmitry Karpeev ystart = 0; 1575af538404SDmitry Karpeev for (j=0; j<Ndomains; ++j) { 1576af538404SDmitry Karpeev maxheight = N/Ndomains + ((N % Ndomains) > j); /* Maximal height of subdomain */ 1577af538404SDmitry 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); 1578af538404SDmitry Karpeev /* Vertical domain limits with an overlap. */ 1579eec7e3faSDmitry Karpeev y[0][0] = PetscMax(ystart - overlap,0); 1580eec7e3faSDmitry Karpeev y[0][1] = PetscMin(ystart + maxheight + overlap,N); 1581eec7e3faSDmitry Karpeev /* Vertical domain limits without an overlap. */ 1582eec7e3faSDmitry Karpeev y[1][0] = ystart; 1583eec7e3faSDmitry Karpeev y[1][1] = PetscMin(ystart + maxheight,N); 1584eec7e3faSDmitry Karpeev xstart = 0; 1585af538404SDmitry Karpeev for (i=0; i<Mdomains; ++i) { 1586af538404SDmitry Karpeev maxwidth = M/Mdomains + ((M % Mdomains) > i); /* Maximal width of subdomain */ 1587af538404SDmitry 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); 1588af538404SDmitry Karpeev /* Horizontal domain limits with an overlap. */ 1589eec7e3faSDmitry Karpeev x[0][0] = PetscMax(xstart - overlap,0); 1590eec7e3faSDmitry Karpeev x[0][1] = PetscMin(xstart + maxwidth + overlap,M); 1591eec7e3faSDmitry Karpeev /* Horizontal domain limits without an overlap. */ 1592eec7e3faSDmitry Karpeev x[1][0] = xstart; 1593eec7e3faSDmitry Karpeev x[1][1] = PetscMin(xstart+maxwidth,M); 15942bbb417fSDmitry Karpeev /* 15952bbb417fSDmitry Karpeev Determine whether this domain intersects this processor's ownership range of pc->pmat. 15962bbb417fSDmitry Karpeev Do this twice: first for the domains with overlaps, and once without. 15972bbb417fSDmitry Karpeev During the first pass create the subcommunicators, and use them on the second pass as well. 15982bbb417fSDmitry Karpeev */ 15992bbb417fSDmitry Karpeev for (q = 0; q < 2; ++q) { 16002bbb417fSDmitry Karpeev /* 16012bbb417fSDmitry Karpeev domain limits, (xleft, xright) and (ylow, yheigh) are adjusted 16022bbb417fSDmitry Karpeev according to whether the domain with an overlap or without is considered. 16032bbb417fSDmitry Karpeev */ 16042bbb417fSDmitry Karpeev xleft = x[q][0]; xright = x[q][1]; 16052bbb417fSDmitry Karpeev ylow = y[q][0]; yhigh = y[q][1]; 1606c3518bceSDmitry Karpeev PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,(&xleft_loc),(&ylow_loc),(&xright_loc),(&yhigh_loc),(&nidx)); 1607af538404SDmitry Karpeev nidx *= dof; 1608eec7e3faSDmitry Karpeev n[q] = nidx; 16092bbb417fSDmitry Karpeev /* 1610eec7e3faSDmitry Karpeev Based on the counted number of indices in the local domain *with an overlap*, 1611af538404SDmitry Karpeev construct a subcommunicator of all the processors supporting this domain. 1612eec7e3faSDmitry Karpeev Observe that a domain with an overlap might have nontrivial local support, 1613eec7e3faSDmitry Karpeev while the domain without an overlap might not. Hence, the decision to participate 1614eec7e3faSDmitry Karpeev in the subcommunicator must be based on the domain with an overlap. 16152bbb417fSDmitry Karpeev */ 16162bbb417fSDmitry Karpeev if (q == 0) { 16172fa5cd67SKarl Rupp if (nidx) color = 1; 16182fa5cd67SKarl Rupp else color = MPI_UNDEFINED; 16192fa5cd67SKarl Rupp 16202bbb417fSDmitry Karpeev ierr = MPI_Comm_split(comm, color, rank, &subcomm);CHKERRQ(ierr); 16212bbb417fSDmitry Karpeev } 1622af538404SDmitry Karpeev /* 1623eec7e3faSDmitry Karpeev Proceed only if the number of local indices *with an overlap* is nonzero. 1624af538404SDmitry Karpeev */ 1625eec7e3faSDmitry Karpeev if (n[0]) { 16262fa5cd67SKarl Rupp if (q == 0) xis = is; 1627af538404SDmitry Karpeev if (q == 1) { 1628af538404SDmitry Karpeev /* 1629eec7e3faSDmitry Karpeev The IS for the no-overlap subdomain shares a communicator with the overlapping domain. 1630af538404SDmitry Karpeev Moreover, if the overlap is zero, the two ISs are identical. 1631af538404SDmitry Karpeev */ 1632b1a0a8a3SJed Brown if (overlap == 0) { 1633eec7e3faSDmitry Karpeev (*is_local)[s] = (*is)[s]; 16342bbb417fSDmitry Karpeev ierr = PetscObjectReference((PetscObject)(*is)[s]);CHKERRQ(ierr); 16352bbb417fSDmitry Karpeev continue; 1636d34fcf5fSBarry Smith } else { 16376a4f0f73SDmitry Karpeev xis = is_local; 1638eec7e3faSDmitry Karpeev subcomm = ((PetscObject)(*is)[s])->comm; 16392bbb417fSDmitry Karpeev } 1640af538404SDmitry Karpeev } /* if (q == 1) */ 16410298fd71SBarry Smith idx = NULL; 1642785e854fSJed Brown ierr = PetscMalloc1(nidx,&idx);CHKERRQ(ierr); 1643eec7e3faSDmitry Karpeev if (nidx) { 16442bbb417fSDmitry Karpeev k = 0; 16452bbb417fSDmitry Karpeev for (jj=ylow_loc; jj<yhigh_loc; ++jj) { 1646af538404SDmitry Karpeev PetscInt x0 = (jj==ylow_loc) ? xleft_loc : xleft; 1647af538404SDmitry Karpeev PetscInt x1 = (jj==yhigh_loc-1) ? xright_loc : xright; 1648af538404SDmitry Karpeev kk = dof*(M*jj + x0); 16492bbb417fSDmitry Karpeev for (ii=x0; ii<x1; ++ii) { 16502bbb417fSDmitry Karpeev for (d = 0; d < dof; ++d) { 16512bbb417fSDmitry Karpeev idx[k++] = kk++; 1652b1a0a8a3SJed Brown } 1653b1a0a8a3SJed Brown } 1654b1a0a8a3SJed Brown } 1655eec7e3faSDmitry Karpeev } 16566a4f0f73SDmitry Karpeev ierr = ISCreateGeneral(subcomm,nidx,idx,PETSC_OWN_POINTER,(*xis)+s);CHKERRQ(ierr); 1657eec7e3faSDmitry Karpeev }/* if (n[0]) */ 16582bbb417fSDmitry Karpeev }/* for (q = 0; q < 2; ++q) */ 16592fa5cd67SKarl Rupp if (n[0]) ++s; 1660af538404SDmitry Karpeev xstart += maxwidth; 1661af538404SDmitry Karpeev } /* for (i = 0; i < Mdomains; ++i) */ 16622bbb417fSDmitry Karpeev ystart += maxheight; 1663af538404SDmitry Karpeev } /* for (j = 0; j < Ndomains; ++j) */ 1664b1a0a8a3SJed Brown PetscFunctionReturn(0); 1665b1a0a8a3SJed Brown } 1666b1a0a8a3SJed Brown 1667b1a0a8a3SJed Brown #undef __FUNCT__ 166806b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMGetSubdomains" 1669b1a0a8a3SJed Brown /*@C 167006b43e7eSDmitry Karpeev PCGASMGetSubdomains - Gets the subdomains supported on this processor 167106b43e7eSDmitry Karpeev for the additive Schwarz preconditioner. 1672b1a0a8a3SJed Brown 1673b1a0a8a3SJed Brown Not Collective 1674b1a0a8a3SJed Brown 1675b1a0a8a3SJed Brown Input Parameter: 1676b1a0a8a3SJed Brown . pc - the preconditioner context 1677b1a0a8a3SJed Brown 1678b1a0a8a3SJed Brown Output Parameters: 1679b1a0a8a3SJed Brown + n - the number of subdomains for this processor (default value = 1) 16800298fd71SBarry Smith . iis - the index sets that define the inner subdomains (without overlap) supported on this processor (can be NULL) 16810298fd71SBarry Smith - ois - the index sets that define the outer subdomains (with overlap) supported on this processor (can be NULL) 1682b1a0a8a3SJed Brown 1683b1a0a8a3SJed Brown 1684b1a0a8a3SJed Brown Notes: 16856a4f0f73SDmitry Karpeev The user is responsible for destroying the ISs and freeing the returned arrays. 1686b1a0a8a3SJed Brown The IS numbering is in the parallel, global numbering of the vector. 1687b1a0a8a3SJed Brown 1688b1a0a8a3SJed Brown Level: advanced 1689b1a0a8a3SJed Brown 169006b43e7eSDmitry Karpeev .keywords: PC, GASM, get, subdomains, additive Schwarz 1691b1a0a8a3SJed Brown 1692*8f3b4b4dSDmitry Karpeev .seealso: PCGASMSetOverlap(), PCGASMGetSubKSP(), PCGASMCreateSubdomains2D(), 1693*8f3b4b4dSDmitry Karpeev PCGASMSetSubdomains(), PCGASMGetSubmatrices() 1694b1a0a8a3SJed Brown @*/ 16956a4f0f73SDmitry Karpeev PetscErrorCode PCGASMGetSubdomains(PC pc,PetscInt *n,IS *iis[],IS *ois[]) 1696b1a0a8a3SJed Brown { 1697f746d493SDmitry Karpeev PC_GASM *osm; 1698b1a0a8a3SJed Brown PetscErrorCode ierr; 1699b1a0a8a3SJed Brown PetscBool match; 17006a4f0f73SDmitry Karpeev PetscInt i; 17015fd66863SKarl Rupp 1702b1a0a8a3SJed Brown PetscFunctionBegin; 1703b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1704251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match);CHKERRQ(ierr); 1705ce94432eSBarry Smith if (!match) SETERRQ2(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Incorrect object type: expected %s, got %s instead", PCGASM, ((PetscObject)pc)->type_name); 1706f746d493SDmitry Karpeev osm = (PC_GASM*)pc->data; 1707ab3e923fSDmitry Karpeev if (n) *n = osm->n; 17086a4f0f73SDmitry Karpeev if (iis) { 1709785e854fSJed Brown ierr = PetscMalloc1(osm->n, iis);CHKERRQ(ierr); 17106a4f0f73SDmitry Karpeev } 17116a4f0f73SDmitry Karpeev if (ois) { 1712785e854fSJed Brown ierr = PetscMalloc1(osm->n, ois);CHKERRQ(ierr); 17136a4f0f73SDmitry Karpeev } 17146a4f0f73SDmitry Karpeev if (iis || ois) { 17156a4f0f73SDmitry Karpeev for (i = 0; i < osm->n; ++i) { 17166a4f0f73SDmitry Karpeev if (iis) (*iis)[i] = osm->iis[i]; 17176a4f0f73SDmitry Karpeev if (ois) (*ois)[i] = osm->ois[i]; 17186a4f0f73SDmitry Karpeev } 1719b1a0a8a3SJed Brown } 1720b1a0a8a3SJed Brown PetscFunctionReturn(0); 1721b1a0a8a3SJed Brown } 1722b1a0a8a3SJed Brown 1723b1a0a8a3SJed Brown #undef __FUNCT__ 172406b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMGetSubmatrices" 1725b1a0a8a3SJed Brown /*@C 172606b43e7eSDmitry Karpeev PCGASMGetSubmatrices - Gets the local submatrices (for this processor 1727b1a0a8a3SJed Brown only) for the additive Schwarz preconditioner. 1728b1a0a8a3SJed Brown 1729b1a0a8a3SJed Brown Not Collective 1730b1a0a8a3SJed Brown 1731b1a0a8a3SJed Brown Input Parameter: 1732b1a0a8a3SJed Brown . pc - the preconditioner context 1733b1a0a8a3SJed Brown 1734b1a0a8a3SJed Brown Output Parameters: 1735b1a0a8a3SJed Brown + n - the number of matrices for this processor (default value = 1) 1736b1a0a8a3SJed Brown - mat - the matrices 1737b1a0a8a3SJed Brown 173806b43e7eSDmitry Karpeev Notes: matrices returned by this routine have the same communicators as the index sets (IS) 1739*8f3b4b4dSDmitry Karpeev used to define subdomains in PCGASMSetSubdomains() 1740b1a0a8a3SJed Brown Level: advanced 1741b1a0a8a3SJed Brown 1742f746d493SDmitry Karpeev .keywords: PC, GASM, set, local, subdomains, additive Schwarz, block Jacobi 1743b1a0a8a3SJed Brown 1744*8f3b4b4dSDmitry Karpeev .seealso: PCGASMSetOverlap(), PCGASMGetSubKSP(), 174506b43e7eSDmitry Karpeev PCGASMCreateSubdomains2D(), PCGASMSetSubdomains(), PCGASMGetSubdomains() 1746b1a0a8a3SJed Brown @*/ 174706b43e7eSDmitry Karpeev PetscErrorCode PCGASMGetSubmatrices(PC pc,PetscInt *n,Mat *mat[]) 1748b1a0a8a3SJed Brown { 1749f746d493SDmitry Karpeev PC_GASM *osm; 1750b1a0a8a3SJed Brown PetscErrorCode ierr; 1751b1a0a8a3SJed Brown PetscBool match; 1752b1a0a8a3SJed Brown 1753b1a0a8a3SJed Brown PetscFunctionBegin; 1754b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1755b1a0a8a3SJed Brown PetscValidIntPointer(n,2); 1756b1a0a8a3SJed Brown if (mat) PetscValidPointer(mat,3); 1757ce94432eSBarry Smith if (!pc->setupcalled) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call after KSPSetUP() or PCSetUp()."); 1758251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match);CHKERRQ(ierr); 1759ce94432eSBarry Smith if (!match) SETERRQ2(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Expected %s, got %s instead", PCGASM, ((PetscObject)pc)->type_name); 1760f746d493SDmitry Karpeev osm = (PC_GASM*)pc->data; 1761ab3e923fSDmitry Karpeev if (n) *n = osm->n; 1762b1a0a8a3SJed Brown if (mat) *mat = osm->pmat; 1763b1a0a8a3SJed Brown PetscFunctionReturn(0); 1764b1a0a8a3SJed Brown } 1765d709ea83SDmitry Karpeev 1766d709ea83SDmitry Karpeev #undef __FUNCT__ 1767*8f3b4b4dSDmitry Karpeev #define __FUNCT__ "PCGASMSetUseDMSubdomains" 1768d709ea83SDmitry Karpeev /*@ 1769*8f3b4b4dSDmitry Karpeev PCGASMSetUseDMSubdomains - Indicates whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible. 1770d709ea83SDmitry Karpeev Logically Collective 1771d709ea83SDmitry Karpeev 1772d709ea83SDmitry Karpeev Input Parameter: 1773d709ea83SDmitry Karpeev + pc - the preconditioner 1774d709ea83SDmitry Karpeev - flg - boolean indicating whether to use subdomains defined by the DM 1775d709ea83SDmitry Karpeev 1776d709ea83SDmitry Karpeev Options Database Key: 1777*8f3b4b4dSDmitry Karpeev . -pc_gasm_dm_subdomains -pc_gasm_overlap -pc_gasm_total_subdomains 1778d709ea83SDmitry Karpeev 1779d709ea83SDmitry Karpeev Level: intermediate 1780d709ea83SDmitry Karpeev 1781d709ea83SDmitry Karpeev Notes: 1782*8f3b4b4dSDmitry Karpeev PCGASMSetSubdomains(), PCGASMSetTotalSubdomains() or PCGASMSetOverlap() take precedence over PCGASMSetUseDMSubdomains(), 1783*8f3b4b4dSDmitry Karpeev so setting PCGASMSetSubdomains() with nontrivial subdomain ISs or any of PCGASMSetTotalSubdomains() and PCGASMSetOverlap() 1784*8f3b4b4dSDmitry Karpeev automatically turns the latter off. 1785d709ea83SDmitry Karpeev 1786d709ea83SDmitry Karpeev .keywords: PC, ASM, DM, set, subdomains, additive Schwarz 1787d709ea83SDmitry Karpeev 1788*8f3b4b4dSDmitry Karpeev .seealso: PCGASMGetUseDMSubdomains(), PCGASMSetSubdomains(), PCGASMSetOverlap() 1789d709ea83SDmitry Karpeev PCGASMCreateSubdomains2D() 1790d709ea83SDmitry Karpeev @*/ 1791*8f3b4b4dSDmitry Karpeev PetscErrorCode PCGASMSetUseDMSubdomains(PC pc,PetscBool flg) 1792d709ea83SDmitry Karpeev { 1793d709ea83SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 1794d709ea83SDmitry Karpeev PetscErrorCode ierr; 1795d709ea83SDmitry Karpeev PetscBool match; 1796d709ea83SDmitry Karpeev 1797d709ea83SDmitry Karpeev PetscFunctionBegin; 1798d709ea83SDmitry Karpeev PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1799d709ea83SDmitry Karpeev PetscValidLogicalCollectiveBool(pc,flg,2); 1800d709ea83SDmitry Karpeev if (pc->setupcalled) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Not for a setup PC."); 1801d709ea83SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match);CHKERRQ(ierr); 1802d709ea83SDmitry Karpeev if (match) { 1803*8f3b4b4dSDmitry Karpeev if (!osm->user_subdomains && osm->N == PETSC_DETERMINE && osm->overlap < 0) { 1804d709ea83SDmitry Karpeev osm->dm_subdomains = flg; 1805d709ea83SDmitry Karpeev } 1806*8f3b4b4dSDmitry Karpeev } 1807d709ea83SDmitry Karpeev PetscFunctionReturn(0); 1808d709ea83SDmitry Karpeev } 1809d709ea83SDmitry Karpeev 1810d709ea83SDmitry Karpeev #undef __FUNCT__ 1811*8f3b4b4dSDmitry Karpeev #define __FUNCT__ "PCGASMGetUseDMSubdomains" 1812d709ea83SDmitry Karpeev /*@ 1813*8f3b4b4dSDmitry Karpeev PCGASMGetUseDMSubdomains - Returns flag indicating whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible. 1814d709ea83SDmitry Karpeev Not Collective 1815d709ea83SDmitry Karpeev 1816d709ea83SDmitry Karpeev Input Parameter: 1817d709ea83SDmitry Karpeev . pc - the preconditioner 1818d709ea83SDmitry Karpeev 1819d709ea83SDmitry Karpeev Output Parameter: 1820d709ea83SDmitry Karpeev . flg - boolean indicating whether to use subdomains defined by the DM 1821d709ea83SDmitry Karpeev 1822d709ea83SDmitry Karpeev Level: intermediate 1823d709ea83SDmitry Karpeev 1824d709ea83SDmitry Karpeev .keywords: PC, ASM, DM, set, subdomains, additive Schwarz 1825d709ea83SDmitry Karpeev 1826*8f3b4b4dSDmitry Karpeev .seealso: PCGASMSetUseDMSubdomains(), PCGASMSetOverlap() 1827d709ea83SDmitry Karpeev PCGASMCreateSubdomains2D() 1828d709ea83SDmitry Karpeev @*/ 1829*8f3b4b4dSDmitry Karpeev PetscErrorCode PCGASMGetUseDMSubdomains(PC pc,PetscBool* flg) 1830d709ea83SDmitry Karpeev { 1831d709ea83SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 1832d709ea83SDmitry Karpeev PetscErrorCode ierr; 1833d709ea83SDmitry Karpeev PetscBool match; 1834d709ea83SDmitry Karpeev 1835d709ea83SDmitry Karpeev PetscFunctionBegin; 1836d709ea83SDmitry Karpeev PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1837d709ea83SDmitry Karpeev PetscValidPointer(flg,2); 1838d709ea83SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match);CHKERRQ(ierr); 1839d709ea83SDmitry Karpeev if (match) { 1840d709ea83SDmitry Karpeev if (flg) *flg = osm->dm_subdomains; 1841d709ea83SDmitry Karpeev } 1842d709ea83SDmitry Karpeev PetscFunctionReturn(0); 1843d709ea83SDmitry Karpeev } 1844