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 76a4f0f73SDmitry Karpeev N - total number of local subdomains on all processors (set in PCGASMSetTotalSubdomains() or calculated in PCSetUp_GASM()) 86a4f0f73SDmitry Karpeev n - actual number of local subdomains on this processor (set in PCGASMSetSubdomains() or calculated in PCGASMSetTotalSubdomains()) 96a4f0f73SDmitry Karpeev nmax - maximum number of local subdomains per processor (calculated in PCGASMSetTotalSubdomains() or in PCSetUp_GASM()) 10b1a0a8a3SJed Brown */ 11b45d2f2cSJed Brown #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 */ 17b1a0a8a3SJed Brown KSP *ksp; /* linear solvers for each block */ 18f746d493SDmitry Karpeev Vec gx,gy; /* Merged work vectors */ 19f746d493SDmitry Karpeev Vec *x,*y; /* Split work vectors; storage aliases pieces of storage of the above merged vectors. */ 206a4f0f73SDmitry Karpeev VecScatter gorestriction; /* merged restriction to disjoint union of outer subdomains */ 216a4f0f73SDmitry Karpeev VecScatter girestriction; /* merged restriction to disjoint union of inner subdomains */ 226a4f0f73SDmitry Karpeev IS *ois; /* index sets that define the outer (conceptually, overlapping) subdomains */ 236a4f0f73SDmitry Karpeev IS *iis; /* index sets that define the inner (conceptually, nonoverlapping) subdomains */ 24f746d493SDmitry Karpeev Mat *pmat; /* subdomain block matrices */ 25f746d493SDmitry Karpeev PCGASMType type; /* use reduced interpolation, restriction or both */ 266a4f0f73SDmitry Karpeev PetscBool create_local; /* whether the autocreated subdomains are local or not. */ 27b1a0a8a3SJed Brown PetscBool type_set; /* if user set this value (so won't change it for symmetric problems) */ 286a4f0f73SDmitry Karpeev PetscBool same_subdomain_solvers; /* flag indicating whether all local solvers are same */ 29b1a0a8a3SJed Brown PetscBool sort_indices; /* flag to sort subdomain indices */ 30d709ea83SDmitry Karpeev PetscBool dm_subdomains; /* whether DM is allowed to define subdomains */ 31f746d493SDmitry Karpeev } PC_GASM; 32b1a0a8a3SJed Brown 33b1a0a8a3SJed Brown #undef __FUNCT__ 3406b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMSubdomainView_Private" 3506b43e7eSDmitry Karpeev static PetscErrorCode PCGASMSubdomainView_Private(PC pc, PetscInt i, PetscViewer viewer) 36af538404SDmitry Karpeev { 37af538404SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 3806b43e7eSDmitry Karpeev PetscInt j,nidx; 39af538404SDmitry Karpeev const PetscInt *idx; 4006b43e7eSDmitry Karpeev PetscViewer sviewer; 41af538404SDmitry Karpeev char *cidx; 42af538404SDmitry Karpeev PetscErrorCode ierr; 43af538404SDmitry Karpeev 44af538404SDmitry Karpeev PetscFunctionBegin; 45ce94432eSBarry 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); 464bde246dSDmitry Karpeev /* Inner subdomains. */ 474bde246dSDmitry Karpeev ierr = ISGetLocalSize(osm->iis[i], &nidx);CHKERRQ(ierr); 484bde246dSDmitry Karpeev /* 494bde246dSDmitry Karpeev No more than 15 characters per index plus a space. 504bde246dSDmitry Karpeev PetscViewerStringSPrintf requires a string of size at least 2, so use (nidx+1) instead of nidx, 514bde246dSDmitry Karpeev in case nidx == 0. That will take care of the space for the trailing '\0' as well. 524bde246dSDmitry Karpeev For nidx == 0, the whole string 16 '\0'. 534bde246dSDmitry Karpeev */ 544bde246dSDmitry Karpeev ierr = PetscMalloc(sizeof(char)*(16*(nidx+1)+1), &cidx);CHKERRQ(ierr); 554bde246dSDmitry Karpeev ierr = PetscViewerStringOpen(PETSC_COMM_SELF, cidx, 16*(nidx+1)+1, &sviewer);CHKERRQ(ierr); 564bde246dSDmitry Karpeev ierr = ISGetIndices(osm->iis[i], &idx);CHKERRQ(ierr); 574bde246dSDmitry Karpeev for (j = 0; j < nidx; ++j) { 584bde246dSDmitry Karpeev ierr = PetscViewerStringSPrintf(sviewer, "%D ", idx[j]);CHKERRQ(ierr); 594bde246dSDmitry Karpeev } 604bde246dSDmitry Karpeev ierr = ISRestoreIndices(osm->iis[i],&idx);CHKERRQ(ierr); 614bde246dSDmitry Karpeev ierr = PetscViewerDestroy(&sviewer);CHKERRQ(ierr); 624bde246dSDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer, "Inner subdomain:\n");CHKERRQ(ierr); 634bde246dSDmitry Karpeev ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 644bde246dSDmitry Karpeev ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr); 654bde246dSDmitry Karpeev ierr = PetscViewerASCIISynchronizedPrintf(viewer, "%s", cidx);CHKERRQ(ierr); 664bde246dSDmitry Karpeev ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 674bde246dSDmitry Karpeev ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr); 684bde246dSDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer, "\n");CHKERRQ(ierr); 694bde246dSDmitry Karpeev ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 704bde246dSDmitry Karpeev ierr = PetscFree(cidx);CHKERRQ(ierr); 714bde246dSDmitry Karpeev /* Outer subdomains. */ 726a4f0f73SDmitry Karpeev ierr = ISGetLocalSize(osm->ois[i], &nidx);CHKERRQ(ierr); 73eec7e3faSDmitry Karpeev /* 74eec7e3faSDmitry Karpeev No more than 15 characters per index plus a space. 75eec7e3faSDmitry Karpeev PetscViewerStringSPrintf requires a string of size at least 2, so use (nidx+1) instead of nidx, 76eec7e3faSDmitry Karpeev in case nidx == 0. That will take care of the space for the trailing '\0' as well. 77eec7e3faSDmitry Karpeev For nidx == 0, the whole string 16 '\0'. 78eec7e3faSDmitry Karpeev */ 79eec7e3faSDmitry Karpeev ierr = PetscMalloc(sizeof(char)*(16*(nidx+1)+1), &cidx);CHKERRQ(ierr); 8006b43e7eSDmitry Karpeev ierr = PetscViewerStringOpen(PETSC_COMM_SELF, cidx, 16*(nidx+1)+1, &sviewer);CHKERRQ(ierr); 816a4f0f73SDmitry Karpeev ierr = ISGetIndices(osm->ois[i], &idx);CHKERRQ(ierr); 82af538404SDmitry Karpeev for (j = 0; j < nidx; ++j) { 83af538404SDmitry Karpeev ierr = PetscViewerStringSPrintf(sviewer,"%D ", idx[j]);CHKERRQ(ierr); 84af538404SDmitry Karpeev } 856bf464f9SBarry Smith ierr = PetscViewerDestroy(&sviewer);CHKERRQ(ierr); 866a4f0f73SDmitry Karpeev ierr = ISRestoreIndices(osm->ois[i],&idx);CHKERRQ(ierr); 876a4f0f73SDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer, "Outer subdomain:\n");CHKERRQ(ierr); 8806b43e7eSDmitry Karpeev ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 897b23a99aSBarry Smith ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr); 90af538404SDmitry Karpeev ierr = PetscViewerASCIISynchronizedPrintf(viewer, "%s", cidx);CHKERRQ(ierr); 91af538404SDmitry Karpeev ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 927b23a99aSBarry Smith ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr); 93af538404SDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer, "\n");CHKERRQ(ierr); 9406b43e7eSDmitry Karpeev ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 95af538404SDmitry Karpeev ierr = PetscFree(cidx);CHKERRQ(ierr); 9606b43e7eSDmitry Karpeev PetscFunctionReturn(0); 9706b43e7eSDmitry Karpeev } 9806b43e7eSDmitry Karpeev 9906b43e7eSDmitry Karpeev #undef __FUNCT__ 10006b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMPrintSubdomains" 10106b43e7eSDmitry Karpeev static PetscErrorCode PCGASMPrintSubdomains(PC pc) 10206b43e7eSDmitry Karpeev { 10306b43e7eSDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 10406b43e7eSDmitry Karpeev const char *prefix; 10506b43e7eSDmitry Karpeev char fname[PETSC_MAX_PATH_LEN+1]; 1064bde246dSDmitry Karpeev PetscInt i, l, d, count, gcount, *permutation, *numbering; 10706b43e7eSDmitry Karpeev PetscBool found; 1080298fd71SBarry Smith PetscViewer viewer, sviewer = NULL; 10906b43e7eSDmitry Karpeev PetscErrorCode ierr; 11006b43e7eSDmitry Karpeev 11106b43e7eSDmitry Karpeev PetscFunctionBegin; 1124bde246dSDmitry Karpeev ierr = PetscMalloc2(osm->n, PetscInt, &permutation, osm->n, PetscInt, &numbering);CHKERRQ(ierr); 1134bde246dSDmitry Karpeev for (i = 0; i < osm->n; ++i) permutation[i] = i; 114ce94432eSBarry Smith ierr = PetscObjectsGetGlobalNumbering(PetscObjectComm((PetscObject)pc), osm->n, (PetscObject*)osm->ois, &gcount, numbering);CHKERRQ(ierr); 1154bde246dSDmitry Karpeev ierr = PetscSortIntWithPermutation(osm->n, numbering, permutation);CHKERRQ(ierr); 11606b43e7eSDmitry Karpeev ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 11706b43e7eSDmitry Karpeev ierr = PetscOptionsGetString(prefix,"-pc_gasm_print_subdomains",fname,PETSC_MAX_PATH_LEN,&found);CHKERRQ(ierr); 11806b43e7eSDmitry Karpeev if (!found) { ierr = PetscStrcpy(fname,"stdout");CHKERRQ(ierr); }; 119ce94432eSBarry Smith ierr = PetscViewerASCIIOpen(PetscObjectComm((PetscObject)pc),fname,&viewer);CHKERRQ(ierr); 1204bde246dSDmitry Karpeev /* 1214bde246dSDmitry Karpeev Make sure the viewer has a name. Otherwise this may cause a deadlock or other weird errors when creating a subcomm viewer: 1224bde246dSDmitry Karpeev the subcomm viewer will attempt to inherit the viewer's name, which, if not set, will be constructed collectively on the comm. 1234bde246dSDmitry Karpeev */ 1244bde246dSDmitry Karpeev ierr = PetscObjectName((PetscObject)viewer);CHKERRQ(ierr); 1254bde246dSDmitry Karpeev l = 0; 1264bde246dSDmitry Karpeev for (count = 0; count < gcount; ++count) { 1274bde246dSDmitry Karpeev /* Now let subdomains go one at a time in the global numbering order and print their subdomain/solver info. */ 1284bde246dSDmitry Karpeev if (l<osm->n) { 1294bde246dSDmitry Karpeev d = permutation[l]; /* d is the local number of the l-th smallest (in the global ordering) among the locally supported subdomains */ 1304bde246dSDmitry Karpeev if (numbering[d] == count) { 1314bde246dSDmitry Karpeev ierr = PetscViewerGetSubcomm(viewer,((PetscObject)osm->ois[d])->comm, &sviewer);CHKERRQ(ierr); 1324bde246dSDmitry Karpeev ierr = PCGASMSubdomainView_Private(pc,d,sviewer);CHKERRQ(ierr); 1334bde246dSDmitry Karpeev ierr = PetscViewerRestoreSubcomm(viewer,((PetscObject)osm->ois[d])->comm, &sviewer);CHKERRQ(ierr); 1344bde246dSDmitry Karpeev ++l; 135af538404SDmitry Karpeev } 1364bde246dSDmitry Karpeev } 137ce94432eSBarry Smith ierr = MPI_Barrier(PetscObjectComm((PetscObject)pc));CHKERRQ(ierr); 1384bde246dSDmitry Karpeev } 1394bde246dSDmitry Karpeev ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 1407105d80bSDmitry Karpeev ierr = PetscFree2(permutation,numbering);CHKERRQ(ierr); 141af538404SDmitry Karpeev PetscFunctionReturn(0); 142af538404SDmitry Karpeev } 143af538404SDmitry Karpeev 144af538404SDmitry Karpeev 145af538404SDmitry Karpeev #undef __FUNCT__ 146f746d493SDmitry Karpeev #define __FUNCT__ "PCView_GASM" 147f746d493SDmitry Karpeev static PetscErrorCode PCView_GASM(PC pc,PetscViewer viewer) 148b1a0a8a3SJed Brown { 149f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 150af538404SDmitry Karpeev const char *prefix; 151b1a0a8a3SJed Brown PetscErrorCode ierr; 152af538404SDmitry Karpeev PetscMPIInt rank, size; 153b1a0a8a3SJed Brown PetscInt i,bsz; 15406b43e7eSDmitry Karpeev PetscBool iascii,view_subdomains=PETSC_FALSE; 155b1a0a8a3SJed Brown PetscViewer sviewer; 15606b43e7eSDmitry Karpeev PetscInt count, l, gcount, *numbering, *permutation; 1576a4f0f73SDmitry Karpeev char overlap[256] = "user-defined overlap"; 1586a4f0f73SDmitry Karpeev char gsubdomains[256] = "unknown total number of subdomains"; 15906b43e7eSDmitry Karpeev char lsubdomains[256] = "unknown number of local subdomains"; 16006b43e7eSDmitry Karpeev char msubdomains[256] = "unknown max number of local subdomains"; 16111aeaf0aSBarry Smith 162b1a0a8a3SJed Brown PetscFunctionBegin; 163ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size);CHKERRQ(ierr); 164ce94432eSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank);CHKERRQ(ierr); 16506b43e7eSDmitry Karpeev ierr = PetscMalloc2(osm->n, PetscInt, &permutation, osm->n, PetscInt, &numbering);CHKERRQ(ierr); 16606b43e7eSDmitry Karpeev for (i = 0; i < osm->n; ++i) permutation[i] = i; 167ce94432eSBarry Smith ierr = PetscObjectsGetGlobalNumbering(PetscObjectComm((PetscObject)pc), osm->n, (PetscObject*)osm->ois, &gcount, numbering);CHKERRQ(ierr); 16806b43e7eSDmitry Karpeev ierr = PetscSortIntWithPermutation(osm->n, numbering, permutation);CHKERRQ(ierr); 16906b43e7eSDmitry Karpeev 17006b43e7eSDmitry Karpeev if (osm->overlap >= 0) { 17106b43e7eSDmitry Karpeev ierr = PetscSNPrintf(overlap,sizeof(overlap),"requested amount of overlap = %D",osm->overlap);CHKERRQ(ierr); 17206b43e7eSDmitry Karpeev } 1736a4f0f73SDmitry Karpeev ierr = PetscSNPrintf(gsubdomains, sizeof(gsubdomains), "total number of subdomains = %D",gcount);CHKERRQ(ierr); 17406b43e7eSDmitry Karpeev if (osm->N > 0) { 17506b43e7eSDmitry Karpeev ierr = PetscSNPrintf(lsubdomains, sizeof(gsubdomains), "number of local subdomains = %D",osm->N);CHKERRQ(ierr); 17606b43e7eSDmitry Karpeev } 17706b43e7eSDmitry Karpeev if (osm->nmax > 0) { 17806b43e7eSDmitry Karpeev ierr = PetscSNPrintf(msubdomains,sizeof(msubdomains),"max number of local subdomains = %D",osm->nmax);CHKERRQ(ierr); 17906b43e7eSDmitry Karpeev } 18006b43e7eSDmitry Karpeev 181af538404SDmitry Karpeev ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 1820298fd71SBarry Smith ierr = PetscOptionsGetBool(prefix,"-pc_gasm_view_subdomains",&view_subdomains,NULL);CHKERRQ(ierr); 18306b43e7eSDmitry Karpeev 18406b43e7eSDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 185b1a0a8a3SJed Brown if (iascii) { 18606b43e7eSDmitry Karpeev /* 18706b43e7eSDmitry Karpeev Make sure the viewer has a name. Otherwise this may cause a deadlock when creating a subcomm viewer: 18806b43e7eSDmitry Karpeev the subcomm viewer will attempt to inherit the viewer's name, which, if not set, will be constructed 18906b43e7eSDmitry Karpeev collectively on the comm. 19006b43e7eSDmitry Karpeev */ 19106b43e7eSDmitry Karpeev ierr = PetscObjectName((PetscObject)viewer);CHKERRQ(ierr); 19206b43e7eSDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer,"Generalized additive Schwarz:\n");CHKERRQ(ierr); 19306b43e7eSDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer,"Restriction/interpolation type: %s\n",PCGASMTypes[osm->type]);CHKERRQ(ierr); 19406b43e7eSDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer,"%s\n",overlap);CHKERRQ(ierr); 19506b43e7eSDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer,"%s\n",gsubdomains);CHKERRQ(ierr); 19606b43e7eSDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer,"%s\n",lsubdomains);CHKERRQ(ierr); 19706b43e7eSDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer,"%s\n",msubdomains);CHKERRQ(ierr); 1987b23a99aSBarry Smith ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr); 199eec7e3faSDmitry Karpeev ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d:%d] number of locally-supported subdomains = %D\n",(int)rank,(int)size,osm->n);CHKERRQ(ierr); 200af538404SDmitry Karpeev ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 2017b23a99aSBarry Smith ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr); 2026a4f0f73SDmitry Karpeev /* Cannot take advantage of osm->same_subdomain_solvers without a global numbering of subdomains. */ 20306b43e7eSDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer,"Subdomain solver info is as follows:\n");CHKERRQ(ierr); 204b1a0a8a3SJed Brown ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 205b1a0a8a3SJed Brown ierr = PetscViewerASCIIPrintf(viewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr); 20606b43e7eSDmitry Karpeev /* Make sure that everybody waits for the banner to be printed. */ 207ce94432eSBarry Smith ierr = MPI_Barrier(PetscObjectComm((PetscObject)viewer));CHKERRQ(ierr); 2084bde246dSDmitry Karpeev /* Now let subdomains go one at a time in the global numbering order and print their subdomain/solver info. */ 20906b43e7eSDmitry Karpeev l = 0; 21006b43e7eSDmitry Karpeev for (count = 0; count < gcount; ++count) { 21106b43e7eSDmitry Karpeev PetscMPIInt srank, ssize; 21206b43e7eSDmitry Karpeev if (l<osm->n) { 21306b43e7eSDmitry Karpeev PetscInt d = permutation[l]; /* d is the local number of the l-th smallest (in the global ordering) among the locally supported subdomains */ 21406b43e7eSDmitry Karpeev if (numbering[d] == count) { 2156a4f0f73SDmitry Karpeev ierr = MPI_Comm_size(((PetscObject)osm->ois[d])->comm, &ssize);CHKERRQ(ierr); 2166a4f0f73SDmitry Karpeev ierr = MPI_Comm_rank(((PetscObject)osm->ois[d])->comm, &srank);CHKERRQ(ierr); 2176a4f0f73SDmitry Karpeev ierr = PetscViewerGetSubcomm(viewer,((PetscObject)osm->ois[d])->comm, &sviewer);CHKERRQ(ierr); 2186a4f0f73SDmitry Karpeev ierr = ISGetLocalSize(osm->ois[d],&bsz);CHKERRQ(ierr); 2197b23a99aSBarry Smith ierr = PetscViewerASCIISynchronizedAllow(sviewer,PETSC_TRUE);CHKERRQ(ierr); 22006b43e7eSDmitry Karpeev ierr = PetscViewerASCIISynchronizedPrintf(sviewer,"[%D:%D] (subcomm [%D:%D]) local subdomain number %D, local size = %D\n",(int)rank,(int)size,(int)srank,(int)ssize,d,bsz);CHKERRQ(ierr); 2216a4f0f73SDmitry Karpeev ierr = PetscViewerFlush(sviewer);CHKERRQ(ierr); 2226a4f0f73SDmitry Karpeev ierr = PetscViewerASCIISynchronizedAllow(sviewer,PETSC_FALSE);CHKERRQ(ierr); 22306b43e7eSDmitry Karpeev if (view_subdomains) { 22406b43e7eSDmitry Karpeev ierr = PCGASMSubdomainView_Private(pc,d,sviewer);CHKERRQ(ierr); 22506b43e7eSDmitry Karpeev } 2266a4f0f73SDmitry Karpeev if (!pc->setupcalled) { 2276a4f0f73SDmitry Karpeev PetscViewerASCIIPrintf(sviewer, "Solver not set up yet: PCSetUp() not yet called\n");CHKERRQ(ierr); 2282fa5cd67SKarl Rupp } else { 22906b43e7eSDmitry Karpeev ierr = KSPView(osm->ksp[d],sviewer);CHKERRQ(ierr); 2306a4f0f73SDmitry Karpeev } 23106b43e7eSDmitry Karpeev ierr = PetscViewerASCIIPrintf(sviewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr); 2327b23a99aSBarry Smith ierr = PetscViewerFlush(sviewer);CHKERRQ(ierr); 2336a4f0f73SDmitry Karpeev ierr = PetscViewerRestoreSubcomm(viewer,((PetscObject)osm->ois[d])->comm, &sviewer);CHKERRQ(ierr); 23406b43e7eSDmitry Karpeev ++l; 235b1a0a8a3SJed Brown } 23606b43e7eSDmitry Karpeev } 237ce94432eSBarry Smith ierr = MPI_Barrier(PetscObjectComm((PetscObject)pc));CHKERRQ(ierr); 238b1a0a8a3SJed Brown } 239b1a0a8a3SJed Brown ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 240b1a0a8a3SJed Brown } 24106b43e7eSDmitry Karpeev ierr = PetscFree2(permutation,numbering);CHKERRQ(ierr); 242b1a0a8a3SJed Brown PetscFunctionReturn(0); 243b1a0a8a3SJed Brown } 244b1a0a8a3SJed Brown 245af538404SDmitry Karpeev 246af538404SDmitry Karpeev 247af538404SDmitry Karpeev 248af538404SDmitry Karpeev 249b1a0a8a3SJed Brown #undef __FUNCT__ 250f746d493SDmitry Karpeev #define __FUNCT__ "PCSetUp_GASM" 251f746d493SDmitry Karpeev static PetscErrorCode PCSetUp_GASM(PC pc) 252b1a0a8a3SJed Brown { 253f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 254b1a0a8a3SJed Brown PetscErrorCode ierr; 255b1a0a8a3SJed Brown PetscBool symset,flg; 25688389755SJed Brown PetscInt i; 257c10bc1a1SDmitry Karpeev PetscMPIInt rank, size; 258b1a0a8a3SJed Brown MatReuse scall = MAT_REUSE_MATRIX; 259b1a0a8a3SJed Brown KSP ksp; 260b1a0a8a3SJed Brown PC subpc; 261b1a0a8a3SJed Brown const char *prefix,*pprefix; 262f746d493SDmitry Karpeev Vec x,y; 2636a4f0f73SDmitry Karpeev PetscInt oni; /* Number of indices in the i-th local outer subdomain. */ 2646a4f0f73SDmitry Karpeev const PetscInt *oidxi; /* Indices from the i-th subdomain local outer subdomain. */ 2656a4f0f73SDmitry Karpeev PetscInt on; /* Number of indices in the disjoint union of local outer subdomains. */ 2666a4f0f73SDmitry Karpeev PetscInt *oidx; /* Indices in the disjoint union of local outer subdomains. */ 2676a4f0f73SDmitry Karpeev IS gois; /* Disjoint union the global indices of outer subdomains. */ 2686a4f0f73SDmitry Karpeev IS goid; /* Identity IS of the size of the disjoint union of outer subdomains. */ 269f746d493SDmitry Karpeev PetscScalar *gxarray, *gyarray; 2702fa5cd67SKarl Rupp PetscInt gofirst; /* Start of locally-owned indices in the vectors -- osm->gx,osm->gy -- over the disjoint union of outer subdomains. */ 2710298fd71SBarry Smith DM *subdomain_dm = NULL; 272b1a0a8a3SJed Brown 273b1a0a8a3SJed Brown PetscFunctionBegin; 274ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr); 275ce94432eSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);CHKERRQ(ierr); 276b1a0a8a3SJed Brown if (!pc->setupcalled) { 277b1a0a8a3SJed Brown 278b1a0a8a3SJed Brown if (!osm->type_set) { 279b1a0a8a3SJed Brown ierr = MatIsSymmetricKnown(pc->pmat,&symset,&flg);CHKERRQ(ierr); 2802fa5cd67SKarl Rupp if (symset && flg) osm->type = PC_GASM_BASIC; 281b1a0a8a3SJed Brown } 282b1a0a8a3SJed Brown 28306b43e7eSDmitry Karpeev /* 28406b43e7eSDmitry Karpeev If subdomains have been set, then the local number of subdomains, osm->n, is NOT PETSC_DECIDE and is at least 1. 28506b43e7eSDmitry Karpeev The total number of subdomains, osm->N is not necessarily set, might be PETSC_DECIDE, and then will have to be calculated from osm->n. 28606b43e7eSDmitry Karpeev */ 28706b43e7eSDmitry Karpeev if (osm->n == PETSC_DECIDE) { 28869ca1f37SDmitry Karpeev /* no subdomains given */ 28965db9045SDmitry Karpeev /* try pc->dm first, if allowed */ 290d709ea83SDmitry Karpeev if (osm->dm_subdomains && pc->dm) { 291a35b7b57SDmitry Karpeev PetscInt num_subdomains, d; 292a35b7b57SDmitry Karpeev char **subdomain_names; 293a35b7b57SDmitry Karpeev IS *inner_subdomain_is, *outer_subdomain_is; 294a35b7b57SDmitry Karpeev ierr = DMCreateDomainDecomposition(pc->dm, &num_subdomains, &subdomain_names, &inner_subdomain_is, &outer_subdomain_is, &subdomain_dm);CHKERRQ(ierr); 295a35b7b57SDmitry Karpeev if (num_subdomains) { 296a35b7b57SDmitry Karpeev ierr = PCGASMSetSubdomains(pc, num_subdomains, inner_subdomain_is, outer_subdomain_is);CHKERRQ(ierr); 29769ca1f37SDmitry Karpeev } 298a35b7b57SDmitry Karpeev for (d = 0; d < num_subdomains; ++d) { 299a35b7b57SDmitry Karpeev if (subdomain_names) {ierr = PetscFree(subdomain_names[d]);CHKERRQ(ierr);} 300a35b7b57SDmitry Karpeev if (inner_subdomain_is) {ierr = ISDestroy(&inner_subdomain_is[d]);CHKERRQ(ierr);} 301a35b7b57SDmitry Karpeev if (outer_subdomain_is) {ierr = ISDestroy(&outer_subdomain_is[d]);CHKERRQ(ierr);} 302b896cd99SPeter Brune if (subdomain_dm) {ierr = DMDestroy(&subdomain_dm[d]);CHKERRQ(ierr);} 303fdc48646SDmitry Karpeev } 304a35b7b57SDmitry Karpeev ierr = PetscFree(subdomain_names);CHKERRQ(ierr); 305a35b7b57SDmitry Karpeev ierr = PetscFree(inner_subdomain_is);CHKERRQ(ierr); 306a35b7b57SDmitry Karpeev ierr = PetscFree(outer_subdomain_is);CHKERRQ(ierr); 307b896cd99SPeter Brune ierr = PetscFree(subdomain_dm);CHKERRQ(ierr); 308fdc48646SDmitry Karpeev } 30906b43e7eSDmitry Karpeev if (osm->n == PETSC_DECIDE) { /* still no subdomains; use one per processor */ 31044fe18e5SDmitry Karpeev osm->nmax = osm->n = 1; 311ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr); 312f746d493SDmitry Karpeev osm->N = size; 313fdc48646SDmitry Karpeev } 31406b43e7eSDmitry Karpeev } 315a35b7b57SDmitry Karpeev if (!osm->iis) { 316a35b7b57SDmitry Karpeev /* 317a35b7b57SDmitry Karpeev The local number of subdomains was set just above, or in PCGASMSetTotalSubdomains(), or in PCGASMSetSubdomains(), 318a35b7b57SDmitry Karpeev but the actual subdomains have not been supplied (in PCGASMSetSubdomains()). 319a35b7b57SDmitry Karpeev We create the requisite number of inner subdomains on PETSC_COMM_SELF (for now). 320a35b7b57SDmitry Karpeev */ 321a35b7b57SDmitry Karpeev ierr = PCGASMCreateLocalSubdomains(pc->pmat,osm->overlap,osm->n,&osm->iis,&osm->ois);CHKERRQ(ierr); 322a35b7b57SDmitry Karpeev } 32306b43e7eSDmitry Karpeev if (osm->N == PETSC_DECIDE) { 3246ac3741eSJed Brown struct {PetscInt max,sum;} inwork,outwork; 325f746d493SDmitry Karpeev /* determine global number of subdomains and the max number of local subdomains */ 3266ac3741eSJed Brown inwork.max = osm->n; 3276ac3741eSJed Brown inwork.sum = osm->n; 328ce94432eSBarry Smith ierr = MPI_Allreduce(&inwork,&outwork,1,MPIU_2INT,PetscMaxSum_Op,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr); 3296ac3741eSJed Brown osm->nmax = outwork.max; 3306ac3741eSJed Brown osm->N = outwork.sum; 331b1a0a8a3SJed Brown } 3326a4f0f73SDmitry Karpeev 333b1a0a8a3SJed Brown ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 334b1a0a8a3SJed Brown flg = PETSC_FALSE; 3350298fd71SBarry Smith ierr = PetscOptionsGetBool(prefix,"-pc_gasm_print_subdomains",&flg,NULL);CHKERRQ(ierr); 336f746d493SDmitry Karpeev if (flg) { ierr = PCGASMPrintSubdomains(pc);CHKERRQ(ierr); } 337b1a0a8a3SJed Brown 338b1a0a8a3SJed Brown if (osm->sort_indices) { 339f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 3406a4f0f73SDmitry Karpeev ierr = ISSort(osm->ois[i]);CHKERRQ(ierr); 3416a4f0f73SDmitry Karpeev ierr = ISSort(osm->iis[i]);CHKERRQ(ierr); 342b1a0a8a3SJed Brown } 343b1a0a8a3SJed Brown } 3446a4f0f73SDmitry Karpeev /* 3456a4f0f73SDmitry Karpeev Merge the ISs, create merged vectors and restrictions. 3466a4f0f73SDmitry Karpeev */ 3476a4f0f73SDmitry Karpeev /* Merge outer subdomain ISs and construct a restriction onto the disjoint union of local outer subdomains. */ 3486a4f0f73SDmitry Karpeev on = 0; 349f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 3506a4f0f73SDmitry Karpeev ierr = ISGetLocalSize(osm->ois[i],&oni);CHKERRQ(ierr); 3516a4f0f73SDmitry Karpeev on += oni; 352f746d493SDmitry Karpeev } 3536a4f0f73SDmitry Karpeev ierr = PetscMalloc(on*sizeof(PetscInt), &oidx);CHKERRQ(ierr); 3546a4f0f73SDmitry Karpeev on = 0; 355f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 3566a4f0f73SDmitry Karpeev ierr = ISGetLocalSize(osm->ois[i],&oni);CHKERRQ(ierr); 3576a4f0f73SDmitry Karpeev ierr = ISGetIndices(osm->ois[i],&oidxi);CHKERRQ(ierr); 3586a4f0f73SDmitry Karpeev ierr = PetscMemcpy(oidx+on, oidxi, sizeof(PetscInt)*oni);CHKERRQ(ierr); 3596a4f0f73SDmitry Karpeev ierr = ISRestoreIndices(osm->ois[i], &oidxi);CHKERRQ(ierr); 3606a4f0f73SDmitry Karpeev on += oni; 361f746d493SDmitry Karpeev } 3626a4f0f73SDmitry Karpeev ierr = ISCreateGeneral(((PetscObject)(pc))->comm, on, oidx, PETSC_OWN_POINTER, &gois);CHKERRQ(ierr); 363f746d493SDmitry Karpeev ierr = MatGetVecs(pc->pmat,&x,&y);CHKERRQ(ierr); 364ce94432eSBarry Smith ierr = VecCreateMPI(PetscObjectComm((PetscObject)pc), on, PETSC_DECIDE, &osm->gx);CHKERRQ(ierr); 365f746d493SDmitry Karpeev ierr = VecDuplicate(osm->gx,&osm->gy);CHKERRQ(ierr); 3660298fd71SBarry Smith ierr = VecGetOwnershipRange(osm->gx, &gofirst, NULL);CHKERRQ(ierr); 367ce94432eSBarry Smith ierr = ISCreateStride(PetscObjectComm((PetscObject)pc),on,gofirst,1, &goid);CHKERRQ(ierr); 3686a4f0f73SDmitry Karpeev ierr = VecScatterCreate(x,gois,osm->gx,goid, &(osm->gorestriction));CHKERRQ(ierr); 3696a4f0f73SDmitry Karpeev ierr = VecDestroy(&x);CHKERRQ(ierr); 3707105d80bSDmitry Karpeev ierr = ISDestroy(&gois);CHKERRQ(ierr); 3712fa5cd67SKarl Rupp 3726a4f0f73SDmitry Karpeev /* Merge inner subdomain ISs and construct a restriction onto the disjoint union of local inner subdomains. */ 3732fa5cd67SKarl Rupp { 3742fa5cd67SKarl Rupp PetscInt ini; /* Number of indices the i-th a local inner subdomain. */ 3756a4f0f73SDmitry Karpeev PetscInt in; /* Number of indices in the disjoint uniont of local inner subdomains. */ 3766a4f0f73SDmitry Karpeev PetscInt *iidx; /* Global indices in the merged local inner subdomain. */ 3776a4f0f73SDmitry Karpeev PetscInt *ioidx; /* Global indices of the disjoint union of inner subdomains within the disjoint union of outer subdomains. */ 3786a4f0f73SDmitry Karpeev IS giis; /* IS for the disjoint union of inner subdomains. */ 3796a4f0f73SDmitry Karpeev IS giois; /* IS for the disjoint union of inner subdomains within the disjoint union of outer subdomains. */ 380f746d493SDmitry Karpeev /**/ 3816a4f0f73SDmitry Karpeev in = 0; 382f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 3836a4f0f73SDmitry Karpeev ierr = ISGetLocalSize(osm->iis[i],&ini);CHKERRQ(ierr); 3846a4f0f73SDmitry Karpeev in += ini; 385f746d493SDmitry Karpeev } 3866a4f0f73SDmitry Karpeev ierr = PetscMalloc(in*sizeof(PetscInt), &iidx);CHKERRQ(ierr); 3876a4f0f73SDmitry Karpeev ierr = PetscMalloc(in*sizeof(PetscInt), &ioidx);CHKERRQ(ierr); 3880298fd71SBarry Smith ierr = VecGetOwnershipRange(osm->gx,&gofirst, NULL);CHKERRQ(ierr); 3896a4f0f73SDmitry Karpeev in = 0; 3906a4f0f73SDmitry Karpeev on = 0; 391f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 3926a4f0f73SDmitry Karpeev const PetscInt *iidxi; /* Global indices of the i-th local inner subdomain. */ 3936a4f0f73SDmitry Karpeev ISLocalToGlobalMapping ltogi; /* Map from global to local indices of the i-th outer local subdomain. */ 3946a4f0f73SDmitry Karpeev PetscInt *ioidxi; /* Local indices of the i-th local inner subdomain within the local outer subdomain. */ 3956a4f0f73SDmitry Karpeev PetscInt ioni; /* Number of indices in ioidxi; if ioni != ini the inner subdomain is not a subdomain of the outer subdomain (error). */ 3966a4f0f73SDmitry Karpeev PetscInt k; 3976a4f0f73SDmitry Karpeev ierr = ISGetLocalSize(osm->iis[i],&ini);CHKERRQ(ierr); 3986a4f0f73SDmitry Karpeev ierr = ISGetLocalSize(osm->ois[i],&oni);CHKERRQ(ierr); 3996a4f0f73SDmitry Karpeev ierr = ISGetIndices(osm->iis[i],&iidxi);CHKERRQ(ierr); 4006a4f0f73SDmitry Karpeev ierr = PetscMemcpy(iidx+in, iidxi, sizeof(PetscInt)*ini);CHKERRQ(ierr); 4016a4f0f73SDmitry Karpeev ierr = ISLocalToGlobalMappingCreateIS(osm->ois[i],<ogi);CHKERRQ(ierr); 4026a4f0f73SDmitry Karpeev ioidxi = ioidx+in; 4036a4f0f73SDmitry Karpeev ierr = ISGlobalToLocalMappingApply(ltogi,IS_GTOLM_DROP,ini,iidxi,&ioni,ioidxi);CHKERRQ(ierr); 4046a4f0f73SDmitry Karpeev ierr = ISLocalToGlobalMappingDestroy(<ogi);CHKERRQ(ierr); 4056a4f0f73SDmitry Karpeev ierr = ISRestoreIndices(osm->iis[i], &iidxi);CHKERRQ(ierr); 4066a4f0f73SDmitry 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); 4072fa5cd67SKarl Rupp for (k = 0; k < ini; ++k) ioidxi[k] += gofirst+on; 4086a4f0f73SDmitry Karpeev in += ini; 4096a4f0f73SDmitry Karpeev on += oni; 410f746d493SDmitry Karpeev } 411ce94432eSBarry Smith ierr = ISCreateGeneral(PetscObjectComm((PetscObject)pc), in, iidx, PETSC_OWN_POINTER, &giis);CHKERRQ(ierr); 412ce94432eSBarry Smith ierr = ISCreateGeneral(PetscObjectComm((PetscObject)pc), in, ioidx, PETSC_OWN_POINTER, &giois);CHKERRQ(ierr); 4136a4f0f73SDmitry Karpeev ierr = VecScatterCreate(y,giis,osm->gy,giois,&osm->girestriction);CHKERRQ(ierr); 4146a4f0f73SDmitry Karpeev ierr = VecDestroy(&y);CHKERRQ(ierr); 4156a4f0f73SDmitry Karpeev ierr = ISDestroy(&giis);CHKERRQ(ierr); 4166a4f0f73SDmitry Karpeev ierr = ISDestroy(&giois);CHKERRQ(ierr); 417b1a0a8a3SJed Brown } 4186a4f0f73SDmitry Karpeev ierr = ISDestroy(&goid);CHKERRQ(ierr); 4192fa5cd67SKarl Rupp 420f746d493SDmitry Karpeev /* Create the subdomain work vectors. */ 421f746d493SDmitry Karpeev ierr = PetscMalloc(osm->n*sizeof(Vec),&osm->x);CHKERRQ(ierr); 422f746d493SDmitry Karpeev ierr = PetscMalloc(osm->n*sizeof(Vec),&osm->y);CHKERRQ(ierr); 423f746d493SDmitry Karpeev ierr = VecGetArray(osm->gx, &gxarray);CHKERRQ(ierr); 424f746d493SDmitry Karpeev ierr = VecGetArray(osm->gy, &gyarray);CHKERRQ(ierr); 4256a4f0f73SDmitry Karpeev for (i=0, on=0; i<osm->n; ++i, on += oni) { 4266a4f0f73SDmitry Karpeev PetscInt oNi; 4276a4f0f73SDmitry Karpeev ierr = ISGetLocalSize(osm->ois[i],&oni);CHKERRQ(ierr); 4286a4f0f73SDmitry Karpeev ierr = ISGetSize(osm->ois[i],&oNi);CHKERRQ(ierr); 4296a4f0f73SDmitry Karpeev ierr = VecCreateMPIWithArray(((PetscObject)(osm->ois[i]))->comm,1,oni,oNi,gxarray+on,&osm->x[i]);CHKERRQ(ierr); 4306a4f0f73SDmitry Karpeev ierr = VecCreateMPIWithArray(((PetscObject)(osm->ois[i]))->comm,1,oni,oNi,gyarray+on,&osm->y[i]);CHKERRQ(ierr); 431b1a0a8a3SJed Brown } 432f746d493SDmitry Karpeev ierr = VecRestoreArray(osm->gx, &gxarray);CHKERRQ(ierr); 433f746d493SDmitry Karpeev ierr = VecRestoreArray(osm->gy, &gyarray);CHKERRQ(ierr); 434b1a0a8a3SJed Brown /* Create the local solvers */ 435f746d493SDmitry Karpeev ierr = PetscMalloc(osm->n*sizeof(KSP*),&osm->ksp);CHKERRQ(ierr); 43644fe18e5SDmitry Karpeev for (i=0; i<osm->n; i++) { /* KSPs are local */ 4376a4f0f73SDmitry Karpeev ierr = KSPCreate(((PetscObject)(osm->ois[i]))->comm,&ksp);CHKERRQ(ierr); 438*3bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ksp);CHKERRQ(ierr); 439b1a0a8a3SJed Brown ierr = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,1);CHKERRQ(ierr); 440b1a0a8a3SJed Brown ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr); 441b1a0a8a3SJed Brown ierr = KSPGetPC(ksp,&subpc);CHKERRQ(ierr); 442b1a0a8a3SJed Brown ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 443b1a0a8a3SJed Brown ierr = KSPSetOptionsPrefix(ksp,prefix);CHKERRQ(ierr); 444b1a0a8a3SJed Brown ierr = KSPAppendOptionsPrefix(ksp,"sub_");CHKERRQ(ierr); 445b1a0a8a3SJed Brown osm->ksp[i] = ksp; 446b1a0a8a3SJed Brown } 447b1a0a8a3SJed Brown scall = MAT_INITIAL_MATRIX; 448b1a0a8a3SJed Brown 4492fa5cd67SKarl Rupp } else { /*if (!pc->setupcalled)*/ 450b1a0a8a3SJed Brown /* 451b1a0a8a3SJed Brown Destroy the blocks from the previous iteration 452b1a0a8a3SJed Brown */ 453b1a0a8a3SJed Brown if (pc->flag == DIFFERENT_NONZERO_PATTERN) { 454f746d493SDmitry Karpeev ierr = MatDestroyMatrices(osm->n,&osm->pmat);CHKERRQ(ierr); 455b1a0a8a3SJed Brown scall = MAT_INITIAL_MATRIX; 456b1a0a8a3SJed Brown } 457b1a0a8a3SJed Brown } 458b1a0a8a3SJed Brown 459b1a0a8a3SJed Brown /* 460f746d493SDmitry Karpeev Extract out the submatrices. 461b1a0a8a3SJed Brown */ 46281a5c4aaSDmitry Karpeev if (size > 1) { 4636a4f0f73SDmitry Karpeev ierr = MatGetSubMatricesParallel(pc->pmat,osm->n,osm->ois, osm->ois,scall,&osm->pmat);CHKERRQ(ierr); 4642fa5cd67SKarl Rupp } else { 4656a4f0f73SDmitry Karpeev ierr = MatGetSubMatrices(pc->pmat,osm->n,osm->ois, osm->ois,scall,&osm->pmat);CHKERRQ(ierr); 46681a5c4aaSDmitry Karpeev } 467b1a0a8a3SJed Brown if (scall == MAT_INITIAL_MATRIX) { 468b1a0a8a3SJed Brown ierr = PetscObjectGetOptionsPrefix((PetscObject)pc->pmat,&pprefix);CHKERRQ(ierr); 469f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 470*3bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)osm->pmat[i]);CHKERRQ(ierr); 471b1a0a8a3SJed Brown ierr = PetscObjectSetOptionsPrefix((PetscObject)osm->pmat[i],pprefix);CHKERRQ(ierr); 472b1a0a8a3SJed Brown } 473b1a0a8a3SJed Brown } 474b1a0a8a3SJed Brown 475b1a0a8a3SJed Brown /* Return control to the user so that the submatrices can be modified (e.g., to apply 476b1a0a8a3SJed Brown different boundary conditions for the submatrices than for the global problem) */ 4776a4f0f73SDmitry Karpeev ierr = PCModifySubMatrices(pc,osm->n,osm->ois,osm->ois,osm->pmat,pc->modifysubmatricesP);CHKERRQ(ierr); 478b1a0a8a3SJed Brown 479b1a0a8a3SJed Brown /* 4806a4f0f73SDmitry Karpeev Loop over submatrices putting them into local ksps 481b1a0a8a3SJed Brown */ 482f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 483b1a0a8a3SJed Brown ierr = KSPSetOperators(osm->ksp[i],osm->pmat[i],osm->pmat[i],pc->flag);CHKERRQ(ierr); 484b1a0a8a3SJed Brown if (!pc->setupcalled) { 485b1a0a8a3SJed Brown ierr = KSPSetFromOptions(osm->ksp[i]);CHKERRQ(ierr); 486b1a0a8a3SJed Brown } 487b1a0a8a3SJed Brown } 488b1a0a8a3SJed Brown PetscFunctionReturn(0); 489b1a0a8a3SJed Brown } 490b1a0a8a3SJed Brown 491b1a0a8a3SJed Brown #undef __FUNCT__ 492f746d493SDmitry Karpeev #define __FUNCT__ "PCSetUpOnBlocks_GASM" 493f746d493SDmitry Karpeev static PetscErrorCode PCSetUpOnBlocks_GASM(PC pc) 494b1a0a8a3SJed Brown { 495f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 496b1a0a8a3SJed Brown PetscErrorCode ierr; 497b1a0a8a3SJed Brown PetscInt i; 498b1a0a8a3SJed Brown 499b1a0a8a3SJed Brown PetscFunctionBegin; 500f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 501b1a0a8a3SJed Brown ierr = KSPSetUp(osm->ksp[i]);CHKERRQ(ierr); 502b1a0a8a3SJed Brown } 503b1a0a8a3SJed Brown PetscFunctionReturn(0); 504b1a0a8a3SJed Brown } 505b1a0a8a3SJed Brown 506b1a0a8a3SJed Brown #undef __FUNCT__ 507f746d493SDmitry Karpeev #define __FUNCT__ "PCApply_GASM" 508f746d493SDmitry Karpeev static PetscErrorCode PCApply_GASM(PC pc,Vec x,Vec y) 509b1a0a8a3SJed Brown { 510f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 511b1a0a8a3SJed Brown PetscErrorCode ierr; 512f746d493SDmitry Karpeev PetscInt i; 513b1a0a8a3SJed Brown ScatterMode forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE; 514b1a0a8a3SJed Brown 515b1a0a8a3SJed Brown PetscFunctionBegin; 516b1a0a8a3SJed Brown /* 5176a4f0f73SDmitry Karpeev Support for limiting the restriction or interpolation only to the inner 518b1a0a8a3SJed Brown subdomain values (leaving the other values 0). 519b1a0a8a3SJed Brown */ 520f746d493SDmitry Karpeev if (!(osm->type & PC_GASM_RESTRICT)) { 521b1a0a8a3SJed Brown /* have to zero the work RHS since scatter may leave some slots empty */ 522f746d493SDmitry Karpeev ierr = VecZeroEntries(osm->gx);CHKERRQ(ierr); 5236a4f0f73SDmitry Karpeev ierr = VecScatterBegin(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 5242fa5cd67SKarl Rupp } else { 5256a4f0f73SDmitry Karpeev ierr = VecScatterBegin(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 526b1a0a8a3SJed Brown } 5276a4f0f73SDmitry Karpeev ierr = VecZeroEntries(osm->gy);CHKERRQ(ierr); 5286a4f0f73SDmitry Karpeev if (!(osm->type & PC_GASM_RESTRICT)) { 5296a4f0f73SDmitry Karpeev ierr = VecScatterEnd(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 5302fa5cd67SKarl Rupp } else { 5316a4f0f73SDmitry Karpeev ierr = VecScatterEnd(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 5326a4f0f73SDmitry Karpeev } 53386b96d36SDmitry Karpeev /* do the subdomain solves */ 53486b96d36SDmitry Karpeev for (i=0; i<osm->n; ++i) { 535b1a0a8a3SJed Brown ierr = KSPSolve(osm->ksp[i],osm->x[i],osm->y[i]);CHKERRQ(ierr); 536b1a0a8a3SJed Brown } 5376a4f0f73SDmitry Karpeev ierr = VecZeroEntries(y);CHKERRQ(ierr); 5386a4f0f73SDmitry Karpeev if (!(osm->type & PC_GASM_INTERPOLATE)) { 5396a4f0f73SDmitry Karpeev ierr = VecScatterBegin(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); 5406a4f0f73SDmitry Karpeev ierr = VecScatterEnd(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); PetscFunctionReturn(0); 5412fa5cd67SKarl Rupp } else { 5426a4f0f73SDmitry Karpeev ierr = VecScatterBegin(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); 5436a4f0f73SDmitry Karpeev ierr = VecScatterEnd(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); PetscFunctionReturn(0); 5446a4f0f73SDmitry Karpeev } 545b1a0a8a3SJed Brown } 546b1a0a8a3SJed Brown 547b1a0a8a3SJed Brown #undef __FUNCT__ 548f746d493SDmitry Karpeev #define __FUNCT__ "PCApplyTranspose_GASM" 549f746d493SDmitry Karpeev static PetscErrorCode PCApplyTranspose_GASM(PC pc,Vec x,Vec y) 550b1a0a8a3SJed Brown { 551f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 552b1a0a8a3SJed Brown PetscErrorCode ierr; 553f746d493SDmitry Karpeev PetscInt i; 554b1a0a8a3SJed Brown ScatterMode forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE; 555b1a0a8a3SJed Brown 556b1a0a8a3SJed Brown PetscFunctionBegin; 557b1a0a8a3SJed Brown /* 558b1a0a8a3SJed Brown Support for limiting the restriction or interpolation to only local 559b1a0a8a3SJed Brown subdomain values (leaving the other values 0). 560b1a0a8a3SJed Brown 561f746d493SDmitry Karpeev Note: these are reversed from the PCApply_GASM() because we are applying the 562b1a0a8a3SJed Brown transpose of the three terms 563b1a0a8a3SJed Brown */ 564f746d493SDmitry Karpeev if (!(osm->type & PC_GASM_INTERPOLATE)) { 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_INTERPOLATE)) { 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 } 577b1a0a8a3SJed Brown /* do the local solves */ 578ab3e923fSDmitry 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. */ 579b1a0a8a3SJed Brown ierr = KSPSolveTranspose(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_RESTRICT)) { 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); 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); 5886a4f0f73SDmitry Karpeev } 589b1a0a8a3SJed Brown PetscFunctionReturn(0); 590b1a0a8a3SJed Brown } 591b1a0a8a3SJed Brown 592b1a0a8a3SJed Brown #undef __FUNCT__ 593a06653b4SBarry Smith #define __FUNCT__ "PCReset_GASM" 594a06653b4SBarry Smith static PetscErrorCode PCReset_GASM(PC pc) 595b1a0a8a3SJed Brown { 596f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 597b1a0a8a3SJed Brown PetscErrorCode ierr; 598b1a0a8a3SJed Brown PetscInt i; 599b1a0a8a3SJed Brown 600b1a0a8a3SJed Brown PetscFunctionBegin; 601b1a0a8a3SJed Brown if (osm->ksp) { 602f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 603a06653b4SBarry Smith ierr = KSPReset(osm->ksp[i]);CHKERRQ(ierr); 604b1a0a8a3SJed Brown } 605b1a0a8a3SJed Brown } 606b1a0a8a3SJed Brown if (osm->pmat) { 607f746d493SDmitry Karpeev if (osm->n > 0) { 608ab3e923fSDmitry Karpeev ierr = MatDestroyMatrices(osm->n,&osm->pmat);CHKERRQ(ierr); 609b1a0a8a3SJed Brown } 610b1a0a8a3SJed Brown } 611d34fcf5fSBarry Smith if (osm->x) { 612f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 613fcfd50ebSBarry Smith ierr = VecDestroy(&osm->x[i]);CHKERRQ(ierr); 614fcfd50ebSBarry Smith ierr = VecDestroy(&osm->y[i]);CHKERRQ(ierr); 615b1a0a8a3SJed Brown } 616d34fcf5fSBarry Smith } 6176bf464f9SBarry Smith ierr = VecDestroy(&osm->gx);CHKERRQ(ierr); 6186bf464f9SBarry Smith ierr = VecDestroy(&osm->gy);CHKERRQ(ierr); 619ab3e923fSDmitry Karpeev 6206a4f0f73SDmitry Karpeev ierr = VecScatterDestroy(&osm->gorestriction);CHKERRQ(ierr); 6216a4f0f73SDmitry Karpeev ierr = VecScatterDestroy(&osm->girestriction);CHKERRQ(ierr); 6226a4f0f73SDmitry Karpeev ierr = PCGASMDestroySubdomains(osm->n,osm->ois,osm->iis);CHKERRQ(ierr); 6236a4f0f73SDmitry Karpeev osm->ois = 0; 6246a4f0f73SDmitry Karpeev osm->iis = 0; 625a06653b4SBarry Smith PetscFunctionReturn(0); 626a06653b4SBarry Smith } 627a06653b4SBarry Smith 628a06653b4SBarry Smith #undef __FUNCT__ 629a06653b4SBarry Smith #define __FUNCT__ "PCDestroy_GASM" 630a06653b4SBarry Smith static PetscErrorCode PCDestroy_GASM(PC pc) 631a06653b4SBarry Smith { 632a06653b4SBarry Smith PC_GASM *osm = (PC_GASM*)pc->data; 633a06653b4SBarry Smith PetscErrorCode ierr; 634a06653b4SBarry Smith PetscInt i; 635a06653b4SBarry Smith 636a06653b4SBarry Smith PetscFunctionBegin; 637a06653b4SBarry Smith ierr = PCReset_GASM(pc);CHKERRQ(ierr); 638a06653b4SBarry Smith if (osm->ksp) { 639a06653b4SBarry Smith for (i=0; i<osm->n; i++) { 6406bf464f9SBarry Smith ierr = KSPDestroy(&osm->ksp[i]);CHKERRQ(ierr); 641a06653b4SBarry Smith } 642a06653b4SBarry Smith ierr = PetscFree(osm->ksp);CHKERRQ(ierr); 643a06653b4SBarry Smith } 644a06653b4SBarry Smith ierr = PetscFree(osm->x);CHKERRQ(ierr); 645a06653b4SBarry Smith ierr = PetscFree(osm->y);CHKERRQ(ierr); 646c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 647b1a0a8a3SJed Brown PetscFunctionReturn(0); 648b1a0a8a3SJed Brown } 649b1a0a8a3SJed Brown 650b1a0a8a3SJed Brown #undef __FUNCT__ 651f746d493SDmitry Karpeev #define __FUNCT__ "PCSetFromOptions_GASM" 652a6dfd86eSKarl Rupp static PetscErrorCode PCSetFromOptions_GASM(PC pc) 653a6dfd86eSKarl Rupp { 654f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 655b1a0a8a3SJed Brown PetscErrorCode ierr; 656b1a0a8a3SJed Brown PetscInt blocks,ovl; 657b1a0a8a3SJed Brown PetscBool symset,flg; 658f746d493SDmitry Karpeev PCGASMType gasmtype; 659b1a0a8a3SJed Brown 660b1a0a8a3SJed Brown PetscFunctionBegin; 661b1a0a8a3SJed Brown /* set the type to symmetric if matrix is symmetric */ 662b1a0a8a3SJed Brown if (!osm->type_set && pc->pmat) { 663b1a0a8a3SJed Brown ierr = MatIsSymmetricKnown(pc->pmat,&symset,&flg);CHKERRQ(ierr); 6642fa5cd67SKarl Rupp if (symset && flg) osm->type = PC_GASM_BASIC; 665b1a0a8a3SJed Brown } 66644fe18e5SDmitry Karpeev ierr = PetscOptionsHead("Generalized additive Schwarz options");CHKERRQ(ierr); 667d709ea83SDmitry Karpeev ierr = PetscOptionsBool("-pc_gasm_dm_subdomains","Use DMCreateDomainDecomposition() to define subdomains","PCGASMSetDMSubdomains",osm->dm_subdomains,&osm->dm_subdomains,&flg);CHKERRQ(ierr); 6686a4f0f73SDmitry Karpeev ierr = PetscOptionsInt("-pc_gasm_total_subdomains","Total number of subdomains across communicator","PCGASMSetTotalSubdomains",osm->n,&blocks,&flg);CHKERRQ(ierr); 66965db9045SDmitry Karpeev if (flg) { 6706a4f0f73SDmitry Karpeev osm->create_local = PETSC_TRUE; 67165db9045SDmitry Karpeev ierr = PetscOptionsBool("-pc_gasm_subdomains_create_local","Whether to make autocreated subdomains local (true by default)","PCGASMSetTotalSubdomains",osm->create_local,&osm->create_local,NULL);CHKERRQ(ierr); 67265db9045SDmitry Karpeev ierr = PCGASMSetTotalSubdomains(pc,blocks,osm->create_local);CHKERRQ(ierr); 673d709ea83SDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 67465db9045SDmitry Karpeev } 67506b43e7eSDmitry Karpeev ierr = PetscOptionsInt("-pc_gasm_overlap","Number of overlapping degrees of freedom","PCGASMSetOverlap",osm->overlap,&ovl,&flg);CHKERRQ(ierr); 67665db9045SDmitry Karpeev if (flg) { 67765db9045SDmitry Karpeev ierr = PCGASMSetOverlap(pc,ovl);CHKERRQ(ierr); 678d709ea83SDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 67965db9045SDmitry Karpeev } 680b1a0a8a3SJed Brown flg = PETSC_FALSE; 681f746d493SDmitry Karpeev ierr = PetscOptionsEnum("-pc_gasm_type","Type of restriction/extension","PCGASMSetType",PCGASMTypes,(PetscEnum)osm->type,(PetscEnum*)&gasmtype,&flg);CHKERRQ(ierr); 682f746d493SDmitry Karpeev if (flg) {ierr = PCGASMSetType(pc,gasmtype);CHKERRQ(ierr); } 683b1a0a8a3SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 684b1a0a8a3SJed Brown PetscFunctionReturn(0); 685b1a0a8a3SJed Brown } 686b1a0a8a3SJed Brown 687b1a0a8a3SJed Brown /*------------------------------------------------------------------------------------*/ 688b1a0a8a3SJed Brown 689b1a0a8a3SJed Brown #undef __FUNCT__ 69006b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMSetSubdomains_GASM" 691f7a08781SBarry Smith static PetscErrorCode PCGASMSetSubdomains_GASM(PC pc,PetscInt n,IS iis[],IS ois[]) 692b1a0a8a3SJed Brown { 693f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 694b1a0a8a3SJed Brown PetscErrorCode ierr; 695b1a0a8a3SJed Brown PetscInt i; 696b1a0a8a3SJed Brown 697b1a0a8a3SJed Brown PetscFunctionBegin; 69869ca1f37SDmitry Karpeev if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Each process must have 1 or more subdomains, n = %D",n); 699ce94432eSBarry Smith if (pc->setupcalled && (n != osm->n || iis || ois)) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetSubdomains() should be called before calling PCSetUp()."); 700b1a0a8a3SJed Brown 701b1a0a8a3SJed Brown if (!pc->setupcalled) { 70293c1ec32SDmitry Karpeev osm->n = n; 7036a4f0f73SDmitry Karpeev osm->ois = 0; 7046a4f0f73SDmitry Karpeev osm->iis = 0; 705a35b7b57SDmitry Karpeev if (ois) { 706a35b7b57SDmitry Karpeev for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)ois[i]);CHKERRQ(ierr);} 707b1a0a8a3SJed Brown } 708a35b7b57SDmitry Karpeev if (iis) { 709a35b7b57SDmitry Karpeev for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)iis[i]);CHKERRQ(ierr);} 710b1a0a8a3SJed Brown } 7116a4f0f73SDmitry Karpeev ierr = PCGASMDestroySubdomains(osm->n,osm->iis,osm->ois);CHKERRQ(ierr); 712a35b7b57SDmitry Karpeev if (ois) { 7136a4f0f73SDmitry Karpeev ierr = PetscMalloc(n*sizeof(IS),&osm->ois);CHKERRQ(ierr); 7142fa5cd67SKarl Rupp for (i=0; i<n; i++) osm->ois[i] = ois[i]; 7156a4f0f73SDmitry Karpeev /* Flag indicating that the user has set outer subdomains, so PCGASM should not increase their size. */ 716b1a0a8a3SJed Brown osm->overlap = -1; 717a35b7b57SDmitry Karpeev if (!iis) { 7186a4f0f73SDmitry Karpeev ierr = PetscMalloc(n*sizeof(IS),&osm->iis);CHKERRQ(ierr); 719a35b7b57SDmitry Karpeev for (i=0; i<n; i++) { 720a35b7b57SDmitry Karpeev for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)ois[i]);CHKERRQ(ierr);} 721a35b7b57SDmitry Karpeev osm->iis[i] = ois[i]; 722a35b7b57SDmitry Karpeev } 723a35b7b57SDmitry Karpeev } 724a35b7b57SDmitry Karpeev } 725a35b7b57SDmitry Karpeev if (iis) { 726a35b7b57SDmitry Karpeev ierr = PetscMalloc(n*sizeof(IS),&osm->iis);CHKERRQ(ierr); 7272fa5cd67SKarl Rupp for (i=0; i<n; i++) osm->iis[i] = iis[i]; 728a35b7b57SDmitry Karpeev if (!ois) { 729a35b7b57SDmitry Karpeev ierr = PetscMalloc(n*sizeof(IS),&osm->ois);CHKERRQ(ierr); 730a35b7b57SDmitry Karpeev for (i=0; i<n; i++) { 731a35b7b57SDmitry Karpeev for (i=0; i<n; i++) { 732a35b7b57SDmitry Karpeev ierr = PetscObjectReference((PetscObject)iis[i]);CHKERRQ(ierr); 733a35b7b57SDmitry Karpeev osm->ois[i] = iis[i]; 734a35b7b57SDmitry Karpeev } 735a35b7b57SDmitry Karpeev } 736a35b7b57SDmitry Karpeev if (osm->overlap > 0) { 737a35b7b57SDmitry Karpeev /* Extend the "overlapping" regions by a number of steps */ 738a35b7b57SDmitry Karpeev ierr = MatIncreaseOverlap(pc->pmat,osm->n,osm->ois,osm->overlap);CHKERRQ(ierr); 739a35b7b57SDmitry Karpeev } 740a35b7b57SDmitry Karpeev } 741b1a0a8a3SJed Brown } 742b1a0a8a3SJed Brown } 743b1a0a8a3SJed Brown PetscFunctionReturn(0); 744b1a0a8a3SJed Brown } 745b1a0a8a3SJed Brown 746b1a0a8a3SJed Brown #undef __FUNCT__ 7476a4f0f73SDmitry Karpeev #define __FUNCT__ "PCGASMSetTotalSubdomains_GASM" 748f7a08781SBarry Smith static PetscErrorCode PCGASMSetTotalSubdomains_GASM(PC pc,PetscInt N, PetscBool create_local) 7490adebc6cSBarry Smith { 750f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 751b1a0a8a3SJed Brown PetscErrorCode ierr; 752b1a0a8a3SJed Brown PetscMPIInt rank,size; 753b1a0a8a3SJed Brown PetscInt n; 7546a4f0f73SDmitry Karpeev PetscInt Nmin, Nmax; 7555fd66863SKarl Rupp 756b1a0a8a3SJed Brown PetscFunctionBegin; 757ce94432eSBarry Smith if (!create_local) SETERRQ(PetscObjectComm((PetscObject)pc), PETSC_ERR_SUP, "No suppor for autocreation of nonlocal subdomains."); 758ce94432eSBarry Smith if (N < 1) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Total number of subdomains must be > 0, N = %D",N); 759ce94432eSBarry Smith ierr = MPI_Allreduce(&N,&Nmin,1,MPIU_INT,MPIU_MIN,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr); 760ce94432eSBarry Smith ierr = MPI_Allreduce(&N,&Nmax,1,MPIU_INT,MPIU_MAX,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr); 761ce94432eSBarry Smith if (Nmin != Nmax) SETERRQ2(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "All processors must use the same number of subdomains. min(N) = %D != %D = max(N)", Nmin, Nmax); 762b1a0a8a3SJed Brown 7636a4f0f73SDmitry Karpeev osm->create_local = create_local; 764b1a0a8a3SJed Brown /* 765b1a0a8a3SJed Brown Split the subdomains equally among all processors 766b1a0a8a3SJed Brown */ 767ce94432eSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);CHKERRQ(ierr); 768ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr); 769b1a0a8a3SJed Brown n = N/size + ((N % size) > rank); 7706a4f0f73SDmitry Karpeev if (!n) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Process %d must have at least one subdomain: total processors %d total blocks %D",(int)rank,(int)size,N); 771f746d493SDmitry Karpeev if (pc->setupcalled && n != osm->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetTotalSubdomains() should be called before PCSetUp()."); 772b1a0a8a3SJed Brown if (!pc->setupcalled) { 7736a4f0f73SDmitry Karpeev ierr = PCGASMDestroySubdomains(osm->n,osm->iis,osm->ois);CHKERRQ(ierr); 77493c1ec32SDmitry Karpeev osm->N = N; 775f746d493SDmitry Karpeev osm->n = n; 77606b43e7eSDmitry Karpeev osm->nmax = N/size + ((N%size) ? 1 : 0); 7776a4f0f73SDmitry Karpeev osm->ois = 0; 7786a4f0f73SDmitry Karpeev osm->iis = 0; 779b1a0a8a3SJed Brown } 780b1a0a8a3SJed Brown PetscFunctionReturn(0); 78106b43e7eSDmitry Karpeev } 782b1a0a8a3SJed Brown 783b1a0a8a3SJed Brown #undef __FUNCT__ 784f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetOverlap_GASM" 785f7a08781SBarry Smith static PetscErrorCode PCGASMSetOverlap_GASM(PC pc,PetscInt ovl) 786b1a0a8a3SJed Brown { 787f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 788b1a0a8a3SJed Brown 789b1a0a8a3SJed Brown PetscFunctionBegin; 790ce94432eSBarry Smith if (ovl < 0) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap value requested"); 791ce94432eSBarry Smith if (pc->setupcalled && ovl != osm->overlap) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetOverlap() should be called before PCSetUp()."); 7922fa5cd67SKarl Rupp if (!pc->setupcalled) osm->overlap = ovl; 793b1a0a8a3SJed Brown PetscFunctionReturn(0); 794b1a0a8a3SJed Brown } 795b1a0a8a3SJed Brown 796b1a0a8a3SJed Brown #undef __FUNCT__ 797f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetType_GASM" 798f7a08781SBarry Smith static PetscErrorCode PCGASMSetType_GASM(PC pc,PCGASMType type) 799b1a0a8a3SJed Brown { 800f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 801b1a0a8a3SJed Brown 802b1a0a8a3SJed Brown PetscFunctionBegin; 803b1a0a8a3SJed Brown osm->type = type; 804b1a0a8a3SJed Brown osm->type_set = PETSC_TRUE; 805b1a0a8a3SJed Brown PetscFunctionReturn(0); 806b1a0a8a3SJed Brown } 807b1a0a8a3SJed Brown 808b1a0a8a3SJed Brown #undef __FUNCT__ 809f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetSortIndices_GASM" 810f7a08781SBarry Smith static PetscErrorCode PCGASMSetSortIndices_GASM(PC pc,PetscBool doSort) 811b1a0a8a3SJed Brown { 812f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 813b1a0a8a3SJed Brown 814b1a0a8a3SJed Brown PetscFunctionBegin; 815b1a0a8a3SJed Brown osm->sort_indices = doSort; 816b1a0a8a3SJed Brown PetscFunctionReturn(0); 817b1a0a8a3SJed Brown } 818b1a0a8a3SJed Brown 819b1a0a8a3SJed Brown #undef __FUNCT__ 820f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMGetSubKSP_GASM" 82144fe18e5SDmitry Karpeev /* 82244fe18e5SDmitry Karpeev FIX: This routine might need to be modified once multiple ranks per subdomain are allowed. 82344fe18e5SDmitry Karpeev In particular, it would upset the global subdomain number calculation. 82444fe18e5SDmitry Karpeev */ 825f7a08781SBarry Smith static PetscErrorCode PCGASMGetSubKSP_GASM(PC pc,PetscInt *n,PetscInt *first,KSP **ksp) 826b1a0a8a3SJed Brown { 827f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 828b1a0a8a3SJed Brown PetscErrorCode ierr; 829b1a0a8a3SJed Brown 830b1a0a8a3SJed Brown PetscFunctionBegin; 831ce94432eSBarry 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"); 832b1a0a8a3SJed Brown 8332fa5cd67SKarl Rupp if (n) *n = osm->n; 834ab3e923fSDmitry Karpeev if (first) { 835ce94432eSBarry Smith ierr = MPI_Scan(&osm->n,first,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr); 836ab3e923fSDmitry Karpeev *first -= osm->n; 837b1a0a8a3SJed Brown } 838b1a0a8a3SJed Brown if (ksp) { 839b1a0a8a3SJed Brown /* Assume that local solves are now different; not necessarily 84006b43e7eSDmitry Karpeev true, though! This flag is used only for PCView_GASM() */ 841b1a0a8a3SJed Brown *ksp = osm->ksp; 8426a4f0f73SDmitry Karpeev osm->same_subdomain_solvers = PETSC_FALSE; 843b1a0a8a3SJed Brown } 844b1a0a8a3SJed Brown PetscFunctionReturn(0); 845ab3e923fSDmitry Karpeev } /* PCGASMGetSubKSP_GASM() */ 846b1a0a8a3SJed Brown 847b1a0a8a3SJed Brown #undef __FUNCT__ 84806b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMSetSubdomains" 849b1a0a8a3SJed Brown /*@C 85006b43e7eSDmitry Karpeev PCGASMSetSubdomains - Sets the subdomains for this processor 85106b43e7eSDmitry Karpeev for the additive Schwarz preconditioner. 852b1a0a8a3SJed Brown 853b1a0a8a3SJed Brown Collective on PC 854b1a0a8a3SJed Brown 855b1a0a8a3SJed Brown Input Parameters: 856b1a0a8a3SJed Brown + pc - the preconditioner context 85706b43e7eSDmitry Karpeev . n - the number of subdomains for this processor 8586a4f0f73SDmitry Karpeev . iis - the index sets that define this processor's local inner subdomains 8590298fd71SBarry Smith (or NULL for PETSc to determine subdomains) 8606a4f0f73SDmitry Karpeev - ois- the index sets that define this processor's local outer subdomains 8610298fd71SBarry Smith (or NULL to use the same as iis) 862b1a0a8a3SJed Brown 863b1a0a8a3SJed Brown Notes: 86406b43e7eSDmitry Karpeev The IS indices use the parallel, global numbering of the vector entries. 8656a4f0f73SDmitry Karpeev Inner subdomains are those where the correction is applied. 8666a4f0f73SDmitry Karpeev Outer subdomains are those where the residual necessary to obtain the 8676a4f0f73SDmitry Karpeev corrections is obtained (see PCGASMType for the use of inner/outer subdomains). 8686a4f0f73SDmitry Karpeev Both inner and outer subdomains can extend over several processors. 8696a4f0f73SDmitry Karpeev This processor's portion of a subdomain is known as a local subdomain. 8706a4f0f73SDmitry Karpeev 8716a4f0f73SDmitry Karpeev By default the GASM preconditioner uses 1 (local) subdomain per processor. 8726a4f0f73SDmitry Karpeev Use PCGASMSetTotalSubdomains() to set the total number of subdomains across 8736a4f0f73SDmitry Karpeev all processors that PCGASM will create automatically, and to specify whether 8746a4f0f73SDmitry Karpeev they should be local or not. 8756a4f0f73SDmitry Karpeev 876b1a0a8a3SJed Brown 877b1a0a8a3SJed Brown Level: advanced 878b1a0a8a3SJed Brown 87906b43e7eSDmitry Karpeev .keywords: PC, GASM, set, subdomains, additive Schwarz 880b1a0a8a3SJed Brown 8816a4f0f73SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetOverlap(), PCGASMGetSubKSP(), 88206b43e7eSDmitry Karpeev PCGASMCreateSubdomains2D(), PCGASMGetSubdomains() 883b1a0a8a3SJed Brown @*/ 8846a4f0f73SDmitry Karpeev PetscErrorCode PCGASMSetSubdomains(PC pc,PetscInt n,IS iis[],IS ois[]) 885b1a0a8a3SJed Brown { 886b1a0a8a3SJed Brown PetscErrorCode ierr; 887b1a0a8a3SJed Brown 888b1a0a8a3SJed Brown PetscFunctionBegin; 889b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 8906a4f0f73SDmitry Karpeev ierr = PetscTryMethod(pc,"PCGASMSetSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,n,iis,ois));CHKERRQ(ierr); 891b1a0a8a3SJed Brown PetscFunctionReturn(0); 892b1a0a8a3SJed Brown } 893b1a0a8a3SJed Brown 894b1a0a8a3SJed Brown #undef __FUNCT__ 8956a4f0f73SDmitry Karpeev #define __FUNCT__ "PCGASMSetTotalSubdomains" 896b1a0a8a3SJed Brown /*@C 8976a4f0f73SDmitry Karpeev PCGASMSetTotalSubdomains - Sets the total number of subdomains to use in the generalized additive 8986a4f0f73SDmitry Karpeev Schwarz preconditioner. The number of subdomains is cumulative across all processors in pc's 8996a4f0f73SDmitry Karpeev communicator. Either all or no processors in the PC communicator must call this routine with 9006a4f0f73SDmitry Karpeev the same N. The subdomains will be created automatically during PCSetUp(). 901b1a0a8a3SJed Brown 902b1a0a8a3SJed Brown Collective on PC 903b1a0a8a3SJed Brown 904b1a0a8a3SJed Brown Input Parameters: 905b1a0a8a3SJed Brown + pc - the preconditioner context 9066a4f0f73SDmitry Karpeev . N - the total number of subdomains cumulative across all processors 9076a4f0f73SDmitry Karpeev - create_local - whether the subdomains to be created are to be local 908b1a0a8a3SJed Brown 909b1a0a8a3SJed Brown Options Database Key: 9106a4f0f73SDmitry Karpeev To set the total number of subdomains and let PCGASM autocreate them, rather than specify the index sets, use the following options: 9116a4f0f73SDmitry Karpeev + -pc_gasm_total_subdomains <n> - sets the total number of subdomains to be autocreated by PCGASM 9126a4f0f73SDmitry Karpeev - -pc_gasm_subdomains_create_local <true|false> - whether autocreated subdomains should be local or not (default is true) 913b1a0a8a3SJed Brown 91406b43e7eSDmitry Karpeev By default the GASM preconditioner uses 1 subdomain per processor. 915b1a0a8a3SJed Brown 9166a4f0f73SDmitry Karpeev 9176a4f0f73SDmitry Karpeev Use PCGASMSetSubdomains() to set subdomains explicitly or to set different numbers 9186a4f0f73SDmitry Karpeev of subdomains per processor. 919b1a0a8a3SJed Brown 920b1a0a8a3SJed Brown Level: advanced 921b1a0a8a3SJed Brown 922f746d493SDmitry Karpeev .keywords: PC, GASM, set, total, global, subdomains, additive Schwarz 923b1a0a8a3SJed Brown 92406b43e7eSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMSetOverlap(), PCGASMGetSubKSP(), 925f746d493SDmitry Karpeev PCGASMCreateSubdomains2D() 926b1a0a8a3SJed Brown @*/ 9276a4f0f73SDmitry Karpeev PetscErrorCode PCGASMSetTotalSubdomains(PC pc,PetscInt N, PetscBool create_local) 928b1a0a8a3SJed Brown { 929b1a0a8a3SJed Brown PetscErrorCode ierr; 930b1a0a8a3SJed Brown 931b1a0a8a3SJed Brown PetscFunctionBegin; 932b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 9336a4f0f73SDmitry Karpeev ierr = PetscTryMethod(pc,"PCGASMSetTotalSubdomains_C",(PC,PetscInt,PetscBool),(pc,N,create_local));CHKERRQ(ierr); 934b1a0a8a3SJed Brown PetscFunctionReturn(0); 935b1a0a8a3SJed Brown } 936b1a0a8a3SJed Brown 937b1a0a8a3SJed Brown #undef __FUNCT__ 938f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetOverlap" 939b1a0a8a3SJed Brown /*@ 940f746d493SDmitry Karpeev PCGASMSetOverlap - Sets the overlap between a pair of subdomains for the 941b1a0a8a3SJed Brown additive Schwarz preconditioner. Either all or no processors in the 942b1a0a8a3SJed Brown PC communicator must call this routine. 943b1a0a8a3SJed Brown 944b1a0a8a3SJed Brown Logically Collective on PC 945b1a0a8a3SJed Brown 946b1a0a8a3SJed Brown Input Parameters: 947b1a0a8a3SJed Brown + pc - the preconditioner context 948b1a0a8a3SJed Brown - ovl - the amount of overlap between subdomains (ovl >= 0, default value = 1) 949b1a0a8a3SJed Brown 950b1a0a8a3SJed Brown Options Database Key: 95106b43e7eSDmitry Karpeev . -pc_gasm_overlap <overlap> - Sets overlap 952b1a0a8a3SJed Brown 953b1a0a8a3SJed Brown Notes: 95406b43e7eSDmitry Karpeev By default the GASM preconditioner uses 1 subdomain per processor. To use 9556a4f0f73SDmitry Karpeev multiple subdomain per perocessor, see PCGASMSetTotalSubdomains() or 9566a4f0f73SDmitry Karpeev PCGASMSetSubdomains() (and the option -pc_gasm_total_subdomains <n>). 957b1a0a8a3SJed Brown 958b1a0a8a3SJed Brown The overlap defaults to 1, so if one desires that no additional 959b1a0a8a3SJed Brown overlap be computed beyond what may have been set with a call to 9606a4f0f73SDmitry Karpeev PCGASMSetTotalSubdomains() or PCGASMSetSubdomains(), then ovl 961b1a0a8a3SJed Brown must be set to be 0. In particular, if one does not explicitly set 96206b43e7eSDmitry Karpeev the subdomains in application code, then all overlap would be computed 963f746d493SDmitry Karpeev internally by PETSc, and using an overlap of 0 would result in an GASM 964b1a0a8a3SJed Brown variant that is equivalent to the block Jacobi preconditioner. 965b1a0a8a3SJed Brown 966b1a0a8a3SJed Brown Note that one can define initial index sets with any overlap via 96706b43e7eSDmitry Karpeev PCGASMSetSubdomains(); the routine PCGASMSetOverlap() merely allows 96806b43e7eSDmitry Karpeev PETSc to extend that overlap further, if desired. 969b1a0a8a3SJed Brown 970b1a0a8a3SJed Brown Level: intermediate 971b1a0a8a3SJed Brown 972f746d493SDmitry Karpeev .keywords: PC, GASM, set, overlap 973b1a0a8a3SJed Brown 9746a4f0f73SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetSubdomains(), PCGASMGetSubKSP(), 97506b43e7eSDmitry Karpeev PCGASMCreateSubdomains2D(), PCGASMGetSubdomains() 976b1a0a8a3SJed Brown @*/ 9777087cfbeSBarry Smith PetscErrorCode PCGASMSetOverlap(PC pc,PetscInt ovl) 978b1a0a8a3SJed Brown { 979b1a0a8a3SJed Brown PetscErrorCode ierr; 980b1a0a8a3SJed Brown 981b1a0a8a3SJed Brown PetscFunctionBegin; 982b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 983b1a0a8a3SJed Brown PetscValidLogicalCollectiveInt(pc,ovl,2); 984f746d493SDmitry Karpeev ierr = PetscTryMethod(pc,"PCGASMSetOverlap_C",(PC,PetscInt),(pc,ovl));CHKERRQ(ierr); 985b1a0a8a3SJed Brown PetscFunctionReturn(0); 986b1a0a8a3SJed Brown } 987b1a0a8a3SJed Brown 988b1a0a8a3SJed Brown #undef __FUNCT__ 989f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetType" 990b1a0a8a3SJed Brown /*@ 991f746d493SDmitry Karpeev PCGASMSetType - Sets the type of restriction and interpolation used 992b1a0a8a3SJed Brown for local problems in the additive Schwarz method. 993b1a0a8a3SJed Brown 994b1a0a8a3SJed Brown Logically Collective on PC 995b1a0a8a3SJed Brown 996b1a0a8a3SJed Brown Input Parameters: 997b1a0a8a3SJed Brown + pc - the preconditioner context 998f746d493SDmitry Karpeev - type - variant of GASM, one of 999b1a0a8a3SJed Brown .vb 1000f746d493SDmitry Karpeev PC_GASM_BASIC - full interpolation and restriction 1001f746d493SDmitry Karpeev PC_GASM_RESTRICT - full restriction, local processor interpolation 1002f746d493SDmitry Karpeev PC_GASM_INTERPOLATE - full interpolation, local processor restriction 1003f746d493SDmitry Karpeev PC_GASM_NONE - local processor restriction and interpolation 1004b1a0a8a3SJed Brown .ve 1005b1a0a8a3SJed Brown 1006b1a0a8a3SJed Brown Options Database Key: 1007f746d493SDmitry Karpeev . -pc_gasm_type [basic,restrict,interpolate,none] - Sets GASM type 1008b1a0a8a3SJed Brown 1009b1a0a8a3SJed Brown Level: intermediate 1010b1a0a8a3SJed Brown 1011f746d493SDmitry Karpeev .keywords: PC, GASM, set, type 1012b1a0a8a3SJed Brown 10136a4f0f73SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetSubdomains(), PCGASMGetSubKSP(), 1014f746d493SDmitry Karpeev PCGASMCreateSubdomains2D() 1015b1a0a8a3SJed Brown @*/ 10167087cfbeSBarry Smith PetscErrorCode PCGASMSetType(PC pc,PCGASMType type) 1017b1a0a8a3SJed Brown { 1018b1a0a8a3SJed Brown PetscErrorCode ierr; 1019b1a0a8a3SJed Brown 1020b1a0a8a3SJed Brown PetscFunctionBegin; 1021b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1022b1a0a8a3SJed Brown PetscValidLogicalCollectiveEnum(pc,type,2); 1023f746d493SDmitry Karpeev ierr = PetscTryMethod(pc,"PCGASMSetType_C",(PC,PCGASMType),(pc,type));CHKERRQ(ierr); 1024b1a0a8a3SJed Brown PetscFunctionReturn(0); 1025b1a0a8a3SJed Brown } 1026b1a0a8a3SJed Brown 1027b1a0a8a3SJed Brown #undef __FUNCT__ 1028f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetSortIndices" 1029b1a0a8a3SJed Brown /*@ 1030f746d493SDmitry Karpeev PCGASMSetSortIndices - Determines whether subdomain indices are sorted. 1031b1a0a8a3SJed Brown 1032b1a0a8a3SJed Brown Logically Collective on PC 1033b1a0a8a3SJed Brown 1034b1a0a8a3SJed Brown Input Parameters: 1035b1a0a8a3SJed Brown + pc - the preconditioner context 1036b1a0a8a3SJed Brown - doSort - sort the subdomain indices 1037b1a0a8a3SJed Brown 1038b1a0a8a3SJed Brown Level: intermediate 1039b1a0a8a3SJed Brown 1040f746d493SDmitry Karpeev .keywords: PC, GASM, set, type 1041b1a0a8a3SJed Brown 10426a4f0f73SDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMSetTotalSubdomains(), PCGASMGetSubKSP(), 1043f746d493SDmitry Karpeev PCGASMCreateSubdomains2D() 1044b1a0a8a3SJed Brown @*/ 10457087cfbeSBarry Smith PetscErrorCode PCGASMSetSortIndices(PC pc,PetscBool doSort) 1046b1a0a8a3SJed Brown { 1047b1a0a8a3SJed Brown PetscErrorCode ierr; 1048b1a0a8a3SJed Brown 1049b1a0a8a3SJed Brown PetscFunctionBegin; 1050b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1051b1a0a8a3SJed Brown PetscValidLogicalCollectiveBool(pc,doSort,2); 1052f746d493SDmitry Karpeev ierr = PetscTryMethod(pc,"PCGASMSetSortIndices_C",(PC,PetscBool),(pc,doSort));CHKERRQ(ierr); 1053b1a0a8a3SJed Brown PetscFunctionReturn(0); 1054b1a0a8a3SJed Brown } 1055b1a0a8a3SJed Brown 1056b1a0a8a3SJed Brown #undef __FUNCT__ 1057f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMGetSubKSP" 1058b1a0a8a3SJed Brown /*@C 1059f746d493SDmitry Karpeev PCGASMGetSubKSP - Gets the local KSP contexts for all blocks on 1060b1a0a8a3SJed Brown this processor. 1061b1a0a8a3SJed Brown 1062b1a0a8a3SJed Brown Collective on PC iff first_local is requested 1063b1a0a8a3SJed Brown 1064b1a0a8a3SJed Brown Input Parameter: 1065b1a0a8a3SJed Brown . pc - the preconditioner context 1066b1a0a8a3SJed Brown 1067b1a0a8a3SJed Brown Output Parameters: 10680298fd71SBarry Smith + n_local - the number of blocks on this processor or NULL 10690298fd71SBarry Smith . first_local - the global number of the first block on this processor or NULL, 10700298fd71SBarry Smith all processors must request or all must pass NULL 1071b1a0a8a3SJed Brown - ksp - the array of KSP contexts 1072b1a0a8a3SJed Brown 1073b1a0a8a3SJed Brown Note: 1074f746d493SDmitry Karpeev After PCGASMGetSubKSP() the array of KSPes is not to be freed 1075b1a0a8a3SJed Brown 1076b1a0a8a3SJed Brown Currently for some matrix implementations only 1 block per processor 1077b1a0a8a3SJed Brown is supported. 1078b1a0a8a3SJed Brown 1079f746d493SDmitry Karpeev You must call KSPSetUp() before calling PCGASMGetSubKSP(). 1080b1a0a8a3SJed Brown 1081b1a0a8a3SJed Brown Level: advanced 1082b1a0a8a3SJed Brown 1083f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, get, sub, KSP, context 1084b1a0a8a3SJed Brown 10856a4f0f73SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetSubdomains(), PCGASMSetOverlap(), 1086f746d493SDmitry Karpeev PCGASMCreateSubdomains2D(), 1087b1a0a8a3SJed Brown @*/ 10887087cfbeSBarry Smith PetscErrorCode PCGASMGetSubKSP(PC pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp[]) 1089b1a0a8a3SJed Brown { 1090b1a0a8a3SJed Brown PetscErrorCode ierr; 1091b1a0a8a3SJed Brown 1092b1a0a8a3SJed Brown PetscFunctionBegin; 1093b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1094f746d493SDmitry Karpeev ierr = PetscUseMethod(pc,"PCGASMGetSubKSP_C",(PC,PetscInt*,PetscInt*,KSP **),(pc,n_local,first_local,ksp));CHKERRQ(ierr); 1095b1a0a8a3SJed Brown PetscFunctionReturn(0); 1096b1a0a8a3SJed Brown } 1097b1a0a8a3SJed Brown 1098b1a0a8a3SJed Brown /* -------------------------------------------------------------------------------------*/ 1099b1a0a8a3SJed Brown /*MC 1100f746d493SDmitry Karpeev PCGASM - Use the (restricted) additive Schwarz method, each block is (approximately) solved with 1101b1a0a8a3SJed Brown its own KSP object. 1102b1a0a8a3SJed Brown 1103b1a0a8a3SJed Brown Options Database Keys: 110406b43e7eSDmitry Karpeev + -pc_gasm_total_block_count <n> - Sets total number of local subdomains (known as blocks) to be distributed among processors 110506b43e7eSDmitry Karpeev . -pc_gasm_view_subdomains - activates the printing of subdomain indices in PCView(), -ksp_view or -snes_view 110606b43e7eSDmitry Karpeev . -pc_gasm_print_subdomains - activates the printing of subdomain indices in PCSetUp() 110706b43e7eSDmitry Karpeev . -pc_gasm_overlap <ovl> - Sets overlap by which to (automatically) extend local subdomains 1108f746d493SDmitry Karpeev - -pc_gasm_type [basic,restrict,interpolate,none] - Sets GASM type 1109b1a0a8a3SJed Brown 1110b1a0a8a3SJed Brown IMPORTANT: If you run with, for example, 3 blocks on 1 processor or 3 blocks on 3 processors you 1111f746d493SDmitry Karpeev will get a different convergence rate due to the default option of -pc_gasm_type restrict. Use 1112f746d493SDmitry Karpeev -pc_gasm_type basic to use the standard GASM. 1113b1a0a8a3SJed Brown 1114b1a0a8a3SJed Brown Notes: Each processor can have one or more blocks, but a block cannot be shared by more 1115b1a0a8a3SJed Brown than one processor. Defaults to one block per processor. 1116b1a0a8a3SJed Brown 1117b1a0a8a3SJed Brown To set options on the solvers for each block append -sub_ to all the KSP, and PC 1118b1a0a8a3SJed Brown options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly 1119b1a0a8a3SJed Brown 1120f746d493SDmitry Karpeev To set the options on the solvers separate for each block call PCGASMGetSubKSP() 1121b1a0a8a3SJed Brown and set the options directly on the resulting KSP object (you can access its PC 1122b1a0a8a3SJed Brown with KSPGetPC()) 1123b1a0a8a3SJed Brown 1124b1a0a8a3SJed Brown 1125b1a0a8a3SJed Brown Level: beginner 1126b1a0a8a3SJed Brown 1127b1a0a8a3SJed Brown Concepts: additive Schwarz method 1128b1a0a8a3SJed Brown 1129b1a0a8a3SJed Brown References: 1130b1a0a8a3SJed Brown An additive variant of the Schwarz alternating method for the case of many subregions 1131b1a0a8a3SJed Brown M Dryja, OB Widlund - Courant Institute, New York University Technical report 1132b1a0a8a3SJed Brown 1133b1a0a8a3SJed Brown Domain Decompositions: Parallel Multilevel Methods for Elliptic Partial Differential Equations, 1134b1a0a8a3SJed Brown Barry Smith, Petter Bjorstad, and William Gropp, Cambridge University Press, ISBN 0-521-49589-X. 1135b1a0a8a3SJed Brown 1136b1a0a8a3SJed Brown .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 113749517cdeSBarry Smith PCBJACOBI, PCGASMGetSubKSP(), PCGASMSetSubdomains(), 11386a4f0f73SDmitry Karpeev PCGASMSetTotalSubdomains(), PCSetModifySubmatrices(), PCGASMSetOverlap(), PCGASMSetType() 1139b1a0a8a3SJed Brown 1140b1a0a8a3SJed Brown M*/ 1141b1a0a8a3SJed Brown 1142b1a0a8a3SJed Brown #undef __FUNCT__ 1143f746d493SDmitry Karpeev #define __FUNCT__ "PCCreate_GASM" 11448cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_GASM(PC pc) 1145b1a0a8a3SJed Brown { 1146b1a0a8a3SJed Brown PetscErrorCode ierr; 1147f746d493SDmitry Karpeev PC_GASM *osm; 1148b1a0a8a3SJed Brown 1149b1a0a8a3SJed Brown PetscFunctionBegin; 1150f746d493SDmitry Karpeev ierr = PetscNewLog(pc,PC_GASM,&osm);CHKERRQ(ierr); 11512fa5cd67SKarl Rupp 1152ab3e923fSDmitry Karpeev osm->N = PETSC_DECIDE; 115306b43e7eSDmitry Karpeev osm->n = PETSC_DECIDE; 1154ab3e923fSDmitry Karpeev osm->nmax = 0; 1155b1a0a8a3SJed Brown osm->overlap = 1; 1156b1a0a8a3SJed Brown osm->ksp = 0; 11576a4f0f73SDmitry Karpeev osm->gorestriction = 0; 11586a4f0f73SDmitry Karpeev osm->girestriction = 0; 1159ab3e923fSDmitry Karpeev osm->gx = 0; 1160ab3e923fSDmitry Karpeev osm->gy = 0; 1161b1a0a8a3SJed Brown osm->x = 0; 1162b1a0a8a3SJed Brown osm->y = 0; 11636a4f0f73SDmitry Karpeev osm->ois = 0; 11646a4f0f73SDmitry Karpeev osm->iis = 0; 1165b1a0a8a3SJed Brown osm->pmat = 0; 1166f746d493SDmitry Karpeev osm->type = PC_GASM_RESTRICT; 11676a4f0f73SDmitry Karpeev osm->same_subdomain_solvers = PETSC_TRUE; 1168b1a0a8a3SJed Brown osm->sort_indices = PETSC_TRUE; 1169d709ea83SDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 1170b1a0a8a3SJed Brown 1171b1a0a8a3SJed Brown pc->data = (void*)osm; 1172f746d493SDmitry Karpeev pc->ops->apply = PCApply_GASM; 1173f746d493SDmitry Karpeev pc->ops->applytranspose = PCApplyTranspose_GASM; 1174f746d493SDmitry Karpeev pc->ops->setup = PCSetUp_GASM; 1175a06653b4SBarry Smith pc->ops->reset = PCReset_GASM; 1176f746d493SDmitry Karpeev pc->ops->destroy = PCDestroy_GASM; 1177f746d493SDmitry Karpeev pc->ops->setfromoptions = PCSetFromOptions_GASM; 1178f746d493SDmitry Karpeev pc->ops->setuponblocks = PCSetUpOnBlocks_GASM; 1179f746d493SDmitry Karpeev pc->ops->view = PCView_GASM; 1180b1a0a8a3SJed Brown pc->ops->applyrichardson = 0; 1181b1a0a8a3SJed Brown 1182bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetSubdomains_C",PCGASMSetSubdomains_GASM);CHKERRQ(ierr); 1183bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetTotalSubdomains_C",PCGASMSetTotalSubdomains_GASM);CHKERRQ(ierr); 1184bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetOverlap_C",PCGASMSetOverlap_GASM);CHKERRQ(ierr); 1185bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetType_C",PCGASMSetType_GASM);CHKERRQ(ierr); 1186bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetSortIndices_C",PCGASMSetSortIndices_GASM);CHKERRQ(ierr); 1187bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGASMGetSubKSP_C",PCGASMGetSubKSP_GASM);CHKERRQ(ierr); 1188b1a0a8a3SJed Brown PetscFunctionReturn(0); 1189b1a0a8a3SJed Brown } 1190b1a0a8a3SJed Brown 1191b1a0a8a3SJed Brown 1192b1a0a8a3SJed Brown #undef __FUNCT__ 119306b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMCreateLocalSubdomains" 1194b1a0a8a3SJed Brown /*@C 119506b43e7eSDmitry Karpeev PCGASMCreateLocalSubdomains - Creates n local index sets for the overlapping 119606b43e7eSDmitry Karpeev Schwarz preconditioner for a any problem based on its matrix. 1197b1a0a8a3SJed Brown 1198b1a0a8a3SJed Brown Collective 1199b1a0a8a3SJed Brown 1200b1a0a8a3SJed Brown Input Parameters: 1201b1a0a8a3SJed Brown + A - The global matrix operator 12026a4f0f73SDmitry Karpeev . overlap - amount of overlap in outer subdomains 120306b43e7eSDmitry Karpeev - n - the number of local subdomains 1204b1a0a8a3SJed Brown 1205b1a0a8a3SJed Brown Output Parameters: 12066a4f0f73SDmitry Karpeev + iis - the array of index sets defining the local inner subdomains (on which the correction is applied) 12076a4f0f73SDmitry Karpeev - ois - the array of index sets defining the local outer subdomains (on which the residual is computed) 1208b1a0a8a3SJed Brown 1209b1a0a8a3SJed Brown Level: advanced 1210b1a0a8a3SJed Brown 12116a4f0f73SDmitry Karpeev Note: this generates n nonoverlapping local inner subdomains on PETSC_COMM_SELF; 12126a4f0f73SDmitry Karpeev PCGASM will generate the overlap from these if you use them in PCGASMSetSubdomains() and set a 12136a4f0f73SDmitry Karpeev nonzero overlap with PCGASMSetOverlap() 1214b1a0a8a3SJed Brown 1215b1a0a8a3SJed Brown In the Fortran version you must provide the array outis[] already allocated of length n. 1216b1a0a8a3SJed Brown 1217f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, create, subdomains, unstructured grid 1218b1a0a8a3SJed Brown 121906b43e7eSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMDestroySubdomains() 1220b1a0a8a3SJed Brown @*/ 12216a4f0f73SDmitry Karpeev PetscErrorCode PCGASMCreateLocalSubdomains(Mat A, PetscInt overlap, PetscInt n, IS *iis[], IS *ois[]) 1222b1a0a8a3SJed Brown { 1223b1a0a8a3SJed Brown MatPartitioning mpart; 1224b1a0a8a3SJed Brown const char *prefix; 122511bd1e4dSLisandro Dalcin PetscErrorCode (*f)(Mat,MatReuse,Mat*); 1226b1a0a8a3SJed Brown PetscMPIInt size; 1227b1a0a8a3SJed Brown PetscInt i,j,rstart,rend,bs; 122811bd1e4dSLisandro Dalcin PetscBool isbaij = PETSC_FALSE,foundpart = PETSC_FALSE; 12290298fd71SBarry Smith Mat Ad = NULL, adj; 1230b1a0a8a3SJed Brown IS ispart,isnumb,*is; 1231b1a0a8a3SJed Brown PetscErrorCode ierr; 1232b1a0a8a3SJed Brown 1233b1a0a8a3SJed Brown PetscFunctionBegin; 1234b1a0a8a3SJed Brown PetscValidHeaderSpecific(A,MAT_CLASSID,1); 12356a4f0f73SDmitry Karpeev PetscValidPointer(iis,4); 1236b1a0a8a3SJed Brown if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of local blocks must be > 0, n = %D",n); 1237b1a0a8a3SJed Brown 1238b1a0a8a3SJed Brown /* Get prefix, row distribution, and block size */ 1239b1a0a8a3SJed Brown ierr = MatGetOptionsPrefix(A,&prefix);CHKERRQ(ierr); 1240b1a0a8a3SJed Brown ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 1241b1a0a8a3SJed Brown ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr); 1242b1a0a8a3SJed 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); 1243b1a0a8a3SJed Brown 1244b1a0a8a3SJed Brown /* Get diagonal block from matrix if possible */ 1245ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr); 12460005d66cSJed Brown ierr = PetscObjectQueryFunction((PetscObject)A,"MatGetDiagonalBlock_C",&f);CHKERRQ(ierr); 1247b1a0a8a3SJed Brown if (f) { 124811bd1e4dSLisandro Dalcin ierr = MatGetDiagonalBlock(A,&Ad);CHKERRQ(ierr); 1249b1a0a8a3SJed Brown } else if (size == 1) { 125011bd1e4dSLisandro Dalcin Ad = A; 1251b1a0a8a3SJed Brown } 1252b1a0a8a3SJed Brown if (Ad) { 1253251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)Ad,MATSEQBAIJ,&isbaij);CHKERRQ(ierr); 1254251f4c67SDmitry Karpeev if (!isbaij) {ierr = PetscObjectTypeCompare((PetscObject)Ad,MATSEQSBAIJ,&isbaij);CHKERRQ(ierr);} 1255b1a0a8a3SJed Brown } 1256b1a0a8a3SJed Brown if (Ad && n > 1) { 1257b1a0a8a3SJed Brown PetscBool match,done; 1258b1a0a8a3SJed Brown /* Try to setup a good matrix partitioning if available */ 1259b1a0a8a3SJed Brown ierr = MatPartitioningCreate(PETSC_COMM_SELF,&mpart);CHKERRQ(ierr); 1260b1a0a8a3SJed Brown ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr); 1261b1a0a8a3SJed Brown ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr); 1262251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGCURRENT,&match);CHKERRQ(ierr); 1263b1a0a8a3SJed Brown if (!match) { 1264251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGSQUARE,&match);CHKERRQ(ierr); 1265b1a0a8a3SJed Brown } 1266b1a0a8a3SJed Brown if (!match) { /* assume a "good" partitioner is available */ 12671a83f524SJed Brown PetscInt na; 12681a83f524SJed Brown const PetscInt *ia,*ja; 1269b1a0a8a3SJed Brown ierr = MatGetRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr); 1270b1a0a8a3SJed Brown if (done) { 1271b1a0a8a3SJed Brown /* Build adjacency matrix by hand. Unfortunately a call to 1272b1a0a8a3SJed Brown MatConvert(Ad,MATMPIADJ,MAT_INITIAL_MATRIX,&adj) will 1273b1a0a8a3SJed Brown remove the block-aij structure and we cannot expect 1274b1a0a8a3SJed Brown MatPartitioning to split vertices as we need */ 12751a83f524SJed Brown PetscInt i,j,len,nnz,cnt,*iia=0,*jja=0; 12761a83f524SJed Brown const PetscInt *row; 1277b1a0a8a3SJed Brown nnz = 0; 1278b1a0a8a3SJed Brown for (i=0; i<na; i++) { /* count number of nonzeros */ 1279b1a0a8a3SJed Brown len = ia[i+1] - ia[i]; 1280b1a0a8a3SJed Brown row = ja + ia[i]; 1281b1a0a8a3SJed Brown for (j=0; j<len; j++) { 1282b1a0a8a3SJed Brown if (row[j] == i) { /* don't count diagonal */ 1283b1a0a8a3SJed Brown len--; break; 1284b1a0a8a3SJed Brown } 1285b1a0a8a3SJed Brown } 1286b1a0a8a3SJed Brown nnz += len; 1287b1a0a8a3SJed Brown } 1288b1a0a8a3SJed Brown ierr = PetscMalloc((na+1)*sizeof(PetscInt),&iia);CHKERRQ(ierr); 1289b1a0a8a3SJed Brown ierr = PetscMalloc((nnz)*sizeof(PetscInt),&jja);CHKERRQ(ierr); 1290b1a0a8a3SJed Brown nnz = 0; 1291b1a0a8a3SJed Brown iia[0] = 0; 1292b1a0a8a3SJed Brown for (i=0; i<na; i++) { /* fill adjacency */ 1293b1a0a8a3SJed Brown cnt = 0; 1294b1a0a8a3SJed Brown len = ia[i+1] - ia[i]; 1295b1a0a8a3SJed Brown row = ja + ia[i]; 1296b1a0a8a3SJed Brown for (j=0; j<len; j++) { 12972fa5cd67SKarl Rupp if (row[j] != i) jja[nnz+cnt++] = row[j]; /* if not diagonal */ 1298b1a0a8a3SJed Brown } 1299b1a0a8a3SJed Brown nnz += cnt; 1300b1a0a8a3SJed Brown iia[i+1] = nnz; 1301b1a0a8a3SJed Brown } 1302b1a0a8a3SJed Brown /* Partitioning of the adjacency matrix */ 13030298fd71SBarry Smith ierr = MatCreateMPIAdj(PETSC_COMM_SELF,na,na,iia,jja,NULL,&adj);CHKERRQ(ierr); 1304b1a0a8a3SJed Brown ierr = MatPartitioningSetAdjacency(mpart,adj);CHKERRQ(ierr); 1305b1a0a8a3SJed Brown ierr = MatPartitioningSetNParts(mpart,n);CHKERRQ(ierr); 1306b1a0a8a3SJed Brown ierr = MatPartitioningApply(mpart,&ispart);CHKERRQ(ierr); 1307b1a0a8a3SJed Brown ierr = ISPartitioningToNumbering(ispart,&isnumb);CHKERRQ(ierr); 13086bf464f9SBarry Smith ierr = MatDestroy(&adj);CHKERRQ(ierr); 1309b1a0a8a3SJed Brown foundpart = PETSC_TRUE; 1310b1a0a8a3SJed Brown } 1311b1a0a8a3SJed Brown ierr = MatRestoreRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr); 1312b1a0a8a3SJed Brown } 13136bf464f9SBarry Smith ierr = MatPartitioningDestroy(&mpart);CHKERRQ(ierr); 1314b1a0a8a3SJed Brown } 1315b1a0a8a3SJed Brown ierr = PetscMalloc(n*sizeof(IS),&is);CHKERRQ(ierr); 1316b1a0a8a3SJed Brown if (!foundpart) { 1317b1a0a8a3SJed Brown 1318b1a0a8a3SJed Brown /* Partitioning by contiguous chunks of rows */ 1319b1a0a8a3SJed Brown 1320b1a0a8a3SJed Brown PetscInt mbs = (rend-rstart)/bs; 1321b1a0a8a3SJed Brown PetscInt start = rstart; 1322b1a0a8a3SJed Brown for (i=0; i<n; i++) { 1323b1a0a8a3SJed Brown PetscInt count = (mbs/n + ((mbs % n) > i)) * bs; 1324b1a0a8a3SJed Brown ierr = ISCreateStride(PETSC_COMM_SELF,count,start,1,&is[i]);CHKERRQ(ierr); 1325b1a0a8a3SJed Brown start += count; 1326b1a0a8a3SJed Brown } 1327b1a0a8a3SJed Brown 1328b1a0a8a3SJed Brown } else { 1329b1a0a8a3SJed Brown 1330b1a0a8a3SJed Brown /* Partitioning by adjacency of diagonal block */ 1331b1a0a8a3SJed Brown 1332b1a0a8a3SJed Brown const PetscInt *numbering; 1333b1a0a8a3SJed Brown PetscInt *count,nidx,*indices,*newidx,start=0; 1334b1a0a8a3SJed Brown /* Get node count in each partition */ 1335b1a0a8a3SJed Brown ierr = PetscMalloc(n*sizeof(PetscInt),&count);CHKERRQ(ierr); 1336b1a0a8a3SJed Brown ierr = ISPartitioningCount(ispart,n,count);CHKERRQ(ierr); 1337b1a0a8a3SJed Brown if (isbaij && bs > 1) { /* adjust for the block-aij case */ 1338b1a0a8a3SJed Brown for (i=0; i<n; i++) count[i] *= bs; 1339b1a0a8a3SJed Brown } 1340b1a0a8a3SJed Brown /* Build indices from node numbering */ 1341b1a0a8a3SJed Brown ierr = ISGetLocalSize(isnumb,&nidx);CHKERRQ(ierr); 1342b1a0a8a3SJed Brown ierr = PetscMalloc(nidx*sizeof(PetscInt),&indices);CHKERRQ(ierr); 1343b1a0a8a3SJed Brown for (i=0; i<nidx; i++) indices[i] = i; /* needs to be initialized */ 1344b1a0a8a3SJed Brown ierr = ISGetIndices(isnumb,&numbering);CHKERRQ(ierr); 1345b1a0a8a3SJed Brown ierr = PetscSortIntWithPermutation(nidx,numbering,indices);CHKERRQ(ierr); 1346b1a0a8a3SJed Brown ierr = ISRestoreIndices(isnumb,&numbering);CHKERRQ(ierr); 1347b1a0a8a3SJed Brown if (isbaij && bs > 1) { /* adjust for the block-aij case */ 1348b1a0a8a3SJed Brown ierr = PetscMalloc(nidx*bs*sizeof(PetscInt),&newidx);CHKERRQ(ierr); 13492fa5cd67SKarl Rupp for (i=0; i<nidx; i++) { 13502fa5cd67SKarl Rupp for (j=0; j<bs; j++) newidx[i*bs+j] = indices[i]*bs + j; 13512fa5cd67SKarl Rupp } 1352b1a0a8a3SJed Brown ierr = PetscFree(indices);CHKERRQ(ierr); 1353b1a0a8a3SJed Brown nidx *= bs; 1354b1a0a8a3SJed Brown indices = newidx; 1355b1a0a8a3SJed Brown } 1356b1a0a8a3SJed Brown /* Shift to get global indices */ 1357b1a0a8a3SJed Brown for (i=0; i<nidx; i++) indices[i] += rstart; 1358b1a0a8a3SJed Brown 1359b1a0a8a3SJed Brown /* Build the index sets for each block */ 1360b1a0a8a3SJed Brown for (i=0; i<n; i++) { 1361b1a0a8a3SJed Brown ierr = ISCreateGeneral(PETSC_COMM_SELF,count[i],&indices[start],PETSC_COPY_VALUES,&is[i]);CHKERRQ(ierr); 1362b1a0a8a3SJed Brown ierr = ISSort(is[i]);CHKERRQ(ierr); 1363b1a0a8a3SJed Brown start += count[i]; 1364b1a0a8a3SJed Brown } 1365b1a0a8a3SJed Brown 1366b1a0a8a3SJed Brown ierr = PetscFree(count); 1367b1a0a8a3SJed Brown ierr = PetscFree(indices); 1368fcfd50ebSBarry Smith ierr = ISDestroy(&isnumb);CHKERRQ(ierr); 1369fcfd50ebSBarry Smith ierr = ISDestroy(&ispart);CHKERRQ(ierr); 1370b1a0a8a3SJed Brown } 13716a4f0f73SDmitry Karpeev *iis = is; 13726a4f0f73SDmitry Karpeev if (!ois) PetscFunctionReturn(0); 13736a4f0f73SDmitry Karpeev /* 13746a4f0f73SDmitry Karpeev Initially make outer subdomains the same as inner subdomains. If nonzero additional overlap 13756a4f0f73SDmitry Karpeev has been requested, copy the inner subdomains over so they can be modified. 13766a4f0f73SDmitry Karpeev */ 13776a4f0f73SDmitry Karpeev ierr = PetscMalloc(n*sizeof(IS),ois);CHKERRQ(ierr); 13786a4f0f73SDmitry Karpeev for (i=0; i<n; ++i) { 13796a4f0f73SDmitry Karpeev if (overlap > 0) { /* With positive overlap, (*iis)[i] will be modified */ 13806a4f0f73SDmitry Karpeev ierr = ISDuplicate((*iis)[i],(*ois)+i);CHKERRQ(ierr); 13816a4f0f73SDmitry Karpeev ierr = ISCopy((*iis)[i],(*ois)[i]);CHKERRQ(ierr); 13826a4f0f73SDmitry Karpeev } else { 13836a4f0f73SDmitry Karpeev ierr = PetscObjectReference((PetscObject)(*iis)[i]);CHKERRQ(ierr); 13846a4f0f73SDmitry Karpeev (*ois)[i] = (*iis)[i]; 13856a4f0f73SDmitry Karpeev } 13866a4f0f73SDmitry Karpeev } 13876a4f0f73SDmitry Karpeev if (overlap > 0) { 13886a4f0f73SDmitry Karpeev /* Extend the "overlapping" regions by a number of steps */ 13896a4f0f73SDmitry Karpeev ierr = MatIncreaseOverlap(A,n,*ois,overlap);CHKERRQ(ierr); 13906a4f0f73SDmitry Karpeev } 1391b1a0a8a3SJed Brown PetscFunctionReturn(0); 1392b1a0a8a3SJed Brown } 1393b1a0a8a3SJed Brown 1394b1a0a8a3SJed Brown #undef __FUNCT__ 1395f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMDestroySubdomains" 1396b1a0a8a3SJed Brown /*@C 1397f746d493SDmitry Karpeev PCGASMDestroySubdomains - Destroys the index sets created with 139806b43e7eSDmitry Karpeev PCGASMCreateLocalSubdomains() or PCGASMCreateSubdomains2D. Should be 139906b43e7eSDmitry Karpeev called after setting subdomains with PCGASMSetSubdomains(). 1400b1a0a8a3SJed Brown 1401b1a0a8a3SJed Brown Collective 1402b1a0a8a3SJed Brown 1403b1a0a8a3SJed Brown Input Parameters: 1404b1a0a8a3SJed Brown + n - the number of index sets 14056a4f0f73SDmitry Karpeev . iis - the array of inner subdomains, 14060298fd71SBarry Smith - ois - the array of outer subdomains, can be NULL 1407b1a0a8a3SJed Brown 14086a4f0f73SDmitry Karpeev Level: intermediate 14096a4f0f73SDmitry Karpeev 14106a4f0f73SDmitry Karpeev Notes: this is merely a convenience subroutine that walks each list, 14116a4f0f73SDmitry Karpeev destroys each IS on the list, and then frees the list. 1412b1a0a8a3SJed Brown 1413f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, create, subdomains, unstructured grid 1414b1a0a8a3SJed Brown 141506b43e7eSDmitry Karpeev .seealso: PCGASMCreateLocalSubdomains(), PCGASMSetSubdomains() 1416b1a0a8a3SJed Brown @*/ 14176a4f0f73SDmitry Karpeev PetscErrorCode PCGASMDestroySubdomains(PetscInt n, IS iis[], IS ois[]) 1418b1a0a8a3SJed Brown { 1419b1a0a8a3SJed Brown PetscInt i; 1420b1a0a8a3SJed Brown PetscErrorCode ierr; 14215fd66863SKarl Rupp 1422b1a0a8a3SJed Brown PetscFunctionBegin; 1423a35b7b57SDmitry Karpeev if (n <= 0) PetscFunctionReturn(0); 1424a35b7b57SDmitry Karpeev if (iis) { 14256a4f0f73SDmitry Karpeev PetscValidPointer(iis,2); 1426a35b7b57SDmitry Karpeev for (i=0; i<n; i++) { 1427a35b7b57SDmitry Karpeev ierr = ISDestroy(&iis[i]);CHKERRQ(ierr); 1428a35b7b57SDmitry Karpeev } 14296a4f0f73SDmitry Karpeev ierr = PetscFree(iis);CHKERRQ(ierr); 1430a35b7b57SDmitry Karpeev } 14316a4f0f73SDmitry Karpeev if (ois) { 1432a35b7b57SDmitry Karpeev for (i=0; i<n; i++) { 1433a35b7b57SDmitry Karpeev ierr = ISDestroy(&ois[i]);CHKERRQ(ierr); 1434a35b7b57SDmitry Karpeev } 14356a4f0f73SDmitry Karpeev ierr = PetscFree(ois);CHKERRQ(ierr); 1436b1a0a8a3SJed Brown } 1437b1a0a8a3SJed Brown PetscFunctionReturn(0); 1438b1a0a8a3SJed Brown } 1439b1a0a8a3SJed Brown 1440af538404SDmitry Karpeev 1441af538404SDmitry Karpeev #define PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,xleft_loc,ylow_loc,xright_loc,yhigh_loc,n) \ 1442af538404SDmitry Karpeev { \ 1443af538404SDmitry Karpeev PetscInt first_row = first/M, last_row = last/M+1; \ 1444af538404SDmitry Karpeev /* \ 1445af538404SDmitry Karpeev Compute ylow_loc and yhigh_loc so that (ylow_loc,xleft) and (yhigh_loc,xright) are the corners \ 1446af538404SDmitry Karpeev of the bounding box of the intersection of the subdomain with the local ownership range (local \ 1447af538404SDmitry Karpeev subdomain). \ 1448af538404SDmitry Karpeev Also compute xleft_loc and xright_loc as the lower and upper bounds on the first and last rows \ 1449af538404SDmitry Karpeev of the intersection. \ 1450af538404SDmitry Karpeev */ \ 1451af538404SDmitry Karpeev /* ylow_loc is the grid row containing the first element of the local sumbdomain */ \ 1452eec7e3faSDmitry Karpeev *ylow_loc = PetscMax(first_row,ylow); \ 1453af538404SDmitry Karpeev /* xleft_loc is the offset of first element of the local subdomain within its grid row (might actually be outside the local subdomain) */ \ 1454eec7e3faSDmitry Karpeev *xleft_loc = *ylow_loc==first_row ? PetscMax(first%M,xleft) : xleft; \ 1455af538404SDmitry Karpeev /* yhigh_loc is the grid row above the last local subdomain element */ \ 1456eec7e3faSDmitry Karpeev *yhigh_loc = PetscMin(last_row,yhigh); \ 1457af538404SDmitry Karpeev /* xright is the offset of the end of the local subdomain within its grid row (might actually be outside the local subdomain) */ \ 1458eec7e3faSDmitry Karpeev *xright_loc = *yhigh_loc==last_row ? PetscMin(xright,last%M) : xright; \ 1459af538404SDmitry Karpeev /* Now compute the size of the local subdomain n. */ \ 1460c3518bceSDmitry Karpeev *n = 0; \ 1461eec7e3faSDmitry Karpeev if (*ylow_loc < *yhigh_loc) { \ 1462af538404SDmitry Karpeev PetscInt width = xright-xleft; \ 1463eec7e3faSDmitry Karpeev *n += width*(*yhigh_loc-*ylow_loc-1); \ 1464eec7e3faSDmitry Karpeev *n += PetscMin(PetscMax(*xright_loc-xleft,0),width); \ 1465eec7e3faSDmitry Karpeev *n -= PetscMin(PetscMax(*xleft_loc-xleft,0), width); \ 1466af538404SDmitry Karpeev } \ 1467af538404SDmitry Karpeev } 1468af538404SDmitry Karpeev 1469c3518bceSDmitry Karpeev 1470eec7e3faSDmitry Karpeev 1471eec7e3faSDmitry Karpeev #undef __FUNCT__ 1472f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMCreateSubdomains2D" 1473b1a0a8a3SJed Brown /*@ 1474f746d493SDmitry Karpeev PCGASMCreateSubdomains2D - Creates the index sets for the overlapping Schwarz 1475b1a0a8a3SJed Brown preconditioner for a two-dimensional problem on a regular grid. 1476b1a0a8a3SJed Brown 1477af538404SDmitry Karpeev Collective 1478b1a0a8a3SJed Brown 1479b1a0a8a3SJed Brown Input Parameters: 148006b43e7eSDmitry Karpeev + M, N - the global number of grid points in the x and y directions 1481af538404SDmitry Karpeev . Mdomains, Ndomains - the global number of subdomains in the x and y directions 1482b1a0a8a3SJed Brown . dof - degrees of freedom per node 1483b1a0a8a3SJed Brown - overlap - overlap in mesh lines 1484b1a0a8a3SJed Brown 1485b1a0a8a3SJed Brown Output Parameters: 1486af538404SDmitry Karpeev + Nsub - the number of local subdomains created 14876a4f0f73SDmitry Karpeev . iis - array of index sets defining inner (nonoverlapping) subdomains 14886a4f0f73SDmitry Karpeev - ois - array of index sets defining outer (overlapping, if overlap > 0) subdomains 1489b1a0a8a3SJed Brown 1490b1a0a8a3SJed Brown 1491b1a0a8a3SJed Brown Level: advanced 1492b1a0a8a3SJed Brown 1493f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, create, subdomains, 2D, regular grid 1494b1a0a8a3SJed Brown 14956a4f0f73SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetSubdomains(), PCGASMGetSubKSP(), 1496f746d493SDmitry Karpeev PCGASMSetOverlap() 1497b1a0a8a3SJed Brown @*/ 14986a4f0f73SDmitry Karpeev PetscErrorCode PCGASMCreateSubdomains2D(PC pc, PetscInt M,PetscInt N,PetscInt Mdomains,PetscInt Ndomains,PetscInt dof,PetscInt overlap, PetscInt *nsub,IS **iis,IS **ois) 1499b1a0a8a3SJed Brown { 1500b1a0a8a3SJed Brown PetscErrorCode ierr; 15012bbb417fSDmitry Karpeev PetscMPIInt size, rank; 15022bbb417fSDmitry Karpeev PetscInt i, j; 15032bbb417fSDmitry Karpeev PetscInt maxheight, maxwidth; 15042bbb417fSDmitry Karpeev PetscInt xstart, xleft, xright, xleft_loc, xright_loc; 15052bbb417fSDmitry Karpeev PetscInt ystart, ylow, yhigh, ylow_loc, yhigh_loc; 1506eec7e3faSDmitry Karpeev PetscInt x[2][2], y[2][2], n[2]; 15072bbb417fSDmitry Karpeev PetscInt first, last; 15082bbb417fSDmitry Karpeev PetscInt nidx, *idx; 15092bbb417fSDmitry Karpeev PetscInt ii,jj,s,q,d; 1510af538404SDmitry Karpeev PetscInt k,kk; 15112bbb417fSDmitry Karpeev PetscMPIInt color; 15122bbb417fSDmitry Karpeev MPI_Comm comm, subcomm; 15136a4f0f73SDmitry Karpeev IS **xis = 0, **is = ois, **is_local = iis; 1514d34fcf5fSBarry Smith 1515b1a0a8a3SJed Brown PetscFunctionBegin; 15162bbb417fSDmitry Karpeev ierr = PetscObjectGetComm((PetscObject)pc, &comm);CHKERRQ(ierr); 15172bbb417fSDmitry Karpeev ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 15182bbb417fSDmitry Karpeev ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 15192bbb417fSDmitry Karpeev ierr = MatGetOwnershipRange(pc->pmat, &first, &last);CHKERRQ(ierr); 1520d34fcf5fSBarry 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) " 15212bbb417fSDmitry Karpeev "does not respect the number of degrees of freedom per grid point %D", first, last, dof); 1522d34fcf5fSBarry Smith 1523af538404SDmitry Karpeev /* Determine the number of domains with nonzero intersections with the local ownership range. */ 15242bbb417fSDmitry Karpeev s = 0; 1525af538404SDmitry Karpeev ystart = 0; 1526af538404SDmitry Karpeev for (j=0; j<Ndomains; ++j) { 1527af538404SDmitry Karpeev maxheight = N/Ndomains + ((N % Ndomains) > j); /* Maximal height of subdomain */ 1528af538404SDmitry 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); 1529eec7e3faSDmitry Karpeev /* Vertical domain limits with an overlap. */ 1530eec7e3faSDmitry Karpeev ylow = PetscMax(ystart - overlap,0); 1531eec7e3faSDmitry Karpeev yhigh = PetscMin(ystart + maxheight + overlap,N); 1532b1a0a8a3SJed Brown xstart = 0; 1533af538404SDmitry Karpeev for (i=0; i<Mdomains; ++i) { 1534af538404SDmitry Karpeev maxwidth = M/Mdomains + ((M % Mdomains) > i); /* Maximal width of subdomain */ 1535af538404SDmitry 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); 1536eec7e3faSDmitry Karpeev /* Horizontal domain limits with an overlap. */ 1537eec7e3faSDmitry Karpeev xleft = PetscMax(xstart - overlap,0); 1538eec7e3faSDmitry Karpeev xright = PetscMin(xstart + maxwidth + overlap,M); 1539af538404SDmitry Karpeev /* 1540af538404SDmitry Karpeev Determine whether this subdomain intersects this processor's ownership range of pc->pmat. 1541af538404SDmitry Karpeev */ 1542c3518bceSDmitry Karpeev PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,(&xleft_loc),(&ylow_loc),(&xright_loc),(&yhigh_loc),(&nidx)); 15432fa5cd67SKarl Rupp if (nidx) ++s; 1544af538404SDmitry Karpeev xstart += maxwidth; 1545af538404SDmitry Karpeev } /* for (i = 0; i < Mdomains; ++i) */ 1546af538404SDmitry Karpeev ystart += maxheight; 1547af538404SDmitry Karpeev } /* for (j = 0; j < Ndomains; ++j) */ 15482fa5cd67SKarl Rupp 1549af538404SDmitry Karpeev /* Now we can allocate the necessary number of ISs. */ 1550af538404SDmitry Karpeev *nsub = s; 1551af538404SDmitry Karpeev ierr = PetscMalloc((*nsub)*sizeof(IS*),is);CHKERRQ(ierr); 1552af538404SDmitry Karpeev ierr = PetscMalloc((*nsub)*sizeof(IS*),is_local);CHKERRQ(ierr); 1553af538404SDmitry Karpeev s = 0; 1554af538404SDmitry Karpeev ystart = 0; 1555af538404SDmitry Karpeev for (j=0; j<Ndomains; ++j) { 1556af538404SDmitry Karpeev maxheight = N/Ndomains + ((N % Ndomains) > j); /* Maximal height of subdomain */ 1557af538404SDmitry 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); 1558af538404SDmitry Karpeev /* Vertical domain limits with an overlap. */ 1559eec7e3faSDmitry Karpeev y[0][0] = PetscMax(ystart - overlap,0); 1560eec7e3faSDmitry Karpeev y[0][1] = PetscMin(ystart + maxheight + overlap,N); 1561eec7e3faSDmitry Karpeev /* Vertical domain limits without an overlap. */ 1562eec7e3faSDmitry Karpeev y[1][0] = ystart; 1563eec7e3faSDmitry Karpeev y[1][1] = PetscMin(ystart + maxheight,N); 1564eec7e3faSDmitry Karpeev xstart = 0; 1565af538404SDmitry Karpeev for (i=0; i<Mdomains; ++i) { 1566af538404SDmitry Karpeev maxwidth = M/Mdomains + ((M % Mdomains) > i); /* Maximal width of subdomain */ 1567af538404SDmitry 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); 1568af538404SDmitry Karpeev /* Horizontal domain limits with an overlap. */ 1569eec7e3faSDmitry Karpeev x[0][0] = PetscMax(xstart - overlap,0); 1570eec7e3faSDmitry Karpeev x[0][1] = PetscMin(xstart + maxwidth + overlap,M); 1571eec7e3faSDmitry Karpeev /* Horizontal domain limits without an overlap. */ 1572eec7e3faSDmitry Karpeev x[1][0] = xstart; 1573eec7e3faSDmitry Karpeev x[1][1] = PetscMin(xstart+maxwidth,M); 15742bbb417fSDmitry Karpeev /* 15752bbb417fSDmitry Karpeev Determine whether this domain intersects this processor's ownership range of pc->pmat. 15762bbb417fSDmitry Karpeev Do this twice: first for the domains with overlaps, and once without. 15772bbb417fSDmitry Karpeev During the first pass create the subcommunicators, and use them on the second pass as well. 15782bbb417fSDmitry Karpeev */ 15792bbb417fSDmitry Karpeev for (q = 0; q < 2; ++q) { 15802bbb417fSDmitry Karpeev /* 15812bbb417fSDmitry Karpeev domain limits, (xleft, xright) and (ylow, yheigh) are adjusted 15822bbb417fSDmitry Karpeev according to whether the domain with an overlap or without is considered. 15832bbb417fSDmitry Karpeev */ 15842bbb417fSDmitry Karpeev xleft = x[q][0]; xright = x[q][1]; 15852bbb417fSDmitry Karpeev ylow = y[q][0]; yhigh = y[q][1]; 1586c3518bceSDmitry Karpeev PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,(&xleft_loc),(&ylow_loc),(&xright_loc),(&yhigh_loc),(&nidx)); 1587af538404SDmitry Karpeev nidx *= dof; 1588eec7e3faSDmitry Karpeev n[q] = nidx; 15892bbb417fSDmitry Karpeev /* 1590eec7e3faSDmitry Karpeev Based on the counted number of indices in the local domain *with an overlap*, 1591af538404SDmitry Karpeev construct a subcommunicator of all the processors supporting this domain. 1592eec7e3faSDmitry Karpeev Observe that a domain with an overlap might have nontrivial local support, 1593eec7e3faSDmitry Karpeev while the domain without an overlap might not. Hence, the decision to participate 1594eec7e3faSDmitry Karpeev in the subcommunicator must be based on the domain with an overlap. 15952bbb417fSDmitry Karpeev */ 15962bbb417fSDmitry Karpeev if (q == 0) { 15972fa5cd67SKarl Rupp if (nidx) color = 1; 15982fa5cd67SKarl Rupp else color = MPI_UNDEFINED; 15992fa5cd67SKarl Rupp 16002bbb417fSDmitry Karpeev ierr = MPI_Comm_split(comm, color, rank, &subcomm);CHKERRQ(ierr); 16012bbb417fSDmitry Karpeev } 1602af538404SDmitry Karpeev /* 1603eec7e3faSDmitry Karpeev Proceed only if the number of local indices *with an overlap* is nonzero. 1604af538404SDmitry Karpeev */ 1605eec7e3faSDmitry Karpeev if (n[0]) { 16062fa5cd67SKarl Rupp if (q == 0) xis = is; 1607af538404SDmitry Karpeev if (q == 1) { 1608af538404SDmitry Karpeev /* 1609eec7e3faSDmitry Karpeev The IS for the no-overlap subdomain shares a communicator with the overlapping domain. 1610af538404SDmitry Karpeev Moreover, if the overlap is zero, the two ISs are identical. 1611af538404SDmitry Karpeev */ 1612b1a0a8a3SJed Brown if (overlap == 0) { 1613eec7e3faSDmitry Karpeev (*is_local)[s] = (*is)[s]; 16142bbb417fSDmitry Karpeev ierr = PetscObjectReference((PetscObject)(*is)[s]);CHKERRQ(ierr); 16152bbb417fSDmitry Karpeev continue; 1616d34fcf5fSBarry Smith } else { 16176a4f0f73SDmitry Karpeev xis = is_local; 1618eec7e3faSDmitry Karpeev subcomm = ((PetscObject)(*is)[s])->comm; 16192bbb417fSDmitry Karpeev } 1620af538404SDmitry Karpeev } /* if (q == 1) */ 16210298fd71SBarry Smith idx = NULL; 16222bbb417fSDmitry Karpeev ierr = PetscMalloc(nidx*sizeof(PetscInt),&idx);CHKERRQ(ierr); 1623eec7e3faSDmitry Karpeev if (nidx) { 16242bbb417fSDmitry Karpeev k = 0; 16252bbb417fSDmitry Karpeev for (jj=ylow_loc; jj<yhigh_loc; ++jj) { 1626af538404SDmitry Karpeev PetscInt x0 = (jj==ylow_loc) ? xleft_loc : xleft; 1627af538404SDmitry Karpeev PetscInt x1 = (jj==yhigh_loc-1) ? xright_loc : xright; 1628af538404SDmitry Karpeev kk = dof*(M*jj + x0); 16292bbb417fSDmitry Karpeev for (ii=x0; ii<x1; ++ii) { 16302bbb417fSDmitry Karpeev for (d = 0; d < dof; ++d) { 16312bbb417fSDmitry Karpeev idx[k++] = kk++; 1632b1a0a8a3SJed Brown } 1633b1a0a8a3SJed Brown } 1634b1a0a8a3SJed Brown } 1635eec7e3faSDmitry Karpeev } 16366a4f0f73SDmitry Karpeev ierr = ISCreateGeneral(subcomm,nidx,idx,PETSC_OWN_POINTER,(*xis)+s);CHKERRQ(ierr); 1637eec7e3faSDmitry Karpeev }/* if (n[0]) */ 16382bbb417fSDmitry Karpeev }/* for (q = 0; q < 2; ++q) */ 16392fa5cd67SKarl Rupp if (n[0]) ++s; 1640af538404SDmitry Karpeev xstart += maxwidth; 1641af538404SDmitry Karpeev } /* for (i = 0; i < Mdomains; ++i) */ 16422bbb417fSDmitry Karpeev ystart += maxheight; 1643af538404SDmitry Karpeev } /* for (j = 0; j < Ndomains; ++j) */ 1644b1a0a8a3SJed Brown PetscFunctionReturn(0); 1645b1a0a8a3SJed Brown } 1646b1a0a8a3SJed Brown 1647b1a0a8a3SJed Brown #undef __FUNCT__ 164806b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMGetSubdomains" 1649b1a0a8a3SJed Brown /*@C 165006b43e7eSDmitry Karpeev PCGASMGetSubdomains - Gets the subdomains supported on this processor 165106b43e7eSDmitry Karpeev for the additive Schwarz preconditioner. 1652b1a0a8a3SJed Brown 1653b1a0a8a3SJed Brown Not Collective 1654b1a0a8a3SJed Brown 1655b1a0a8a3SJed Brown Input Parameter: 1656b1a0a8a3SJed Brown . pc - the preconditioner context 1657b1a0a8a3SJed Brown 1658b1a0a8a3SJed Brown Output Parameters: 1659b1a0a8a3SJed Brown + n - the number of subdomains for this processor (default value = 1) 16600298fd71SBarry Smith . iis - the index sets that define the inner subdomains (without overlap) supported on this processor (can be NULL) 16610298fd71SBarry Smith - ois - the index sets that define the outer subdomains (with overlap) supported on this processor (can be NULL) 1662b1a0a8a3SJed Brown 1663b1a0a8a3SJed Brown 1664b1a0a8a3SJed Brown Notes: 16656a4f0f73SDmitry Karpeev The user is responsible for destroying the ISs and freeing the returned arrays. 1666b1a0a8a3SJed Brown The IS numbering is in the parallel, global numbering of the vector. 1667b1a0a8a3SJed Brown 1668b1a0a8a3SJed Brown Level: advanced 1669b1a0a8a3SJed Brown 167006b43e7eSDmitry Karpeev .keywords: PC, GASM, get, subdomains, additive Schwarz 1671b1a0a8a3SJed Brown 16726a4f0f73SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetOverlap(), PCGASMGetSubKSP(), 167306b43e7eSDmitry Karpeev PCGASMCreateSubdomains2D(), PCGASMSetSubdomains(), PCGASMGetSubmatrices() 1674b1a0a8a3SJed Brown @*/ 16756a4f0f73SDmitry Karpeev PetscErrorCode PCGASMGetSubdomains(PC pc,PetscInt *n,IS *iis[],IS *ois[]) 1676b1a0a8a3SJed Brown { 1677f746d493SDmitry Karpeev PC_GASM *osm; 1678b1a0a8a3SJed Brown PetscErrorCode ierr; 1679b1a0a8a3SJed Brown PetscBool match; 16806a4f0f73SDmitry Karpeev PetscInt i; 16815fd66863SKarl Rupp 1682b1a0a8a3SJed Brown PetscFunctionBegin; 1683b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1684251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match);CHKERRQ(ierr); 1685ce94432eSBarry Smith if (!match) SETERRQ2(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Incorrect object type: expected %s, got %s instead", PCGASM, ((PetscObject)pc)->type_name); 1686f746d493SDmitry Karpeev osm = (PC_GASM*)pc->data; 1687ab3e923fSDmitry Karpeev if (n) *n = osm->n; 16886a4f0f73SDmitry Karpeev if (iis) { 16896a4f0f73SDmitry Karpeev ierr = PetscMalloc(osm->n*sizeof(IS), iis);CHKERRQ(ierr); 16906a4f0f73SDmitry Karpeev } 16916a4f0f73SDmitry Karpeev if (ois) { 16926a4f0f73SDmitry Karpeev ierr = PetscMalloc(osm->n*sizeof(IS), ois);CHKERRQ(ierr); 16936a4f0f73SDmitry Karpeev } 16946a4f0f73SDmitry Karpeev if (iis || ois) { 16956a4f0f73SDmitry Karpeev for (i = 0; i < osm->n; ++i) { 16966a4f0f73SDmitry Karpeev if (iis) (*iis)[i] = osm->iis[i]; 16976a4f0f73SDmitry Karpeev if (ois) (*ois)[i] = osm->ois[i]; 16986a4f0f73SDmitry Karpeev } 1699b1a0a8a3SJed Brown } 1700b1a0a8a3SJed Brown PetscFunctionReturn(0); 1701b1a0a8a3SJed Brown } 1702b1a0a8a3SJed Brown 1703b1a0a8a3SJed Brown #undef __FUNCT__ 170406b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMGetSubmatrices" 1705b1a0a8a3SJed Brown /*@C 170606b43e7eSDmitry Karpeev PCGASMGetSubmatrices - Gets the local submatrices (for this processor 1707b1a0a8a3SJed Brown only) for the additive Schwarz preconditioner. 1708b1a0a8a3SJed Brown 1709b1a0a8a3SJed Brown Not Collective 1710b1a0a8a3SJed Brown 1711b1a0a8a3SJed Brown Input Parameter: 1712b1a0a8a3SJed Brown . pc - the preconditioner context 1713b1a0a8a3SJed Brown 1714b1a0a8a3SJed Brown Output Parameters: 1715b1a0a8a3SJed Brown + n - the number of matrices for this processor (default value = 1) 1716b1a0a8a3SJed Brown - mat - the matrices 1717b1a0a8a3SJed Brown 171806b43e7eSDmitry Karpeev Notes: matrices returned by this routine have the same communicators as the index sets (IS) 171906b43e7eSDmitry Karpeev used to define subdomains in PCGASMSetSubdomains(), or PETSC_COMM_SELF, if the 17206a4f0f73SDmitry Karpeev subdomains were defined using PCGASMSetTotalSubdomains(). 1721b1a0a8a3SJed Brown Level: advanced 1722b1a0a8a3SJed Brown 1723f746d493SDmitry Karpeev .keywords: PC, GASM, set, local, subdomains, additive Schwarz, block Jacobi 1724b1a0a8a3SJed Brown 17256a4f0f73SDmitry Karpeev .seealso: PCGASMSetTotalSubdomain(), PCGASMSetOverlap(), PCGASMGetSubKSP(), 172606b43e7eSDmitry Karpeev PCGASMCreateSubdomains2D(), PCGASMSetSubdomains(), PCGASMGetSubdomains() 1727b1a0a8a3SJed Brown @*/ 172806b43e7eSDmitry Karpeev PetscErrorCode PCGASMGetSubmatrices(PC pc,PetscInt *n,Mat *mat[]) 1729b1a0a8a3SJed Brown { 1730f746d493SDmitry Karpeev PC_GASM *osm; 1731b1a0a8a3SJed Brown PetscErrorCode ierr; 1732b1a0a8a3SJed Brown PetscBool match; 1733b1a0a8a3SJed Brown 1734b1a0a8a3SJed Brown PetscFunctionBegin; 1735b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1736b1a0a8a3SJed Brown PetscValidIntPointer(n,2); 1737b1a0a8a3SJed Brown if (mat) PetscValidPointer(mat,3); 1738ce94432eSBarry Smith if (!pc->setupcalled) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call after KSPSetUP() or PCSetUp()."); 1739251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match);CHKERRQ(ierr); 1740ce94432eSBarry Smith if (!match) SETERRQ2(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Expected %s, got %s instead", PCGASM, ((PetscObject)pc)->type_name); 1741f746d493SDmitry Karpeev osm = (PC_GASM*)pc->data; 1742ab3e923fSDmitry Karpeev if (n) *n = osm->n; 1743b1a0a8a3SJed Brown if (mat) *mat = osm->pmat; 1744b1a0a8a3SJed Brown PetscFunctionReturn(0); 1745b1a0a8a3SJed Brown } 1746d709ea83SDmitry Karpeev 1747d709ea83SDmitry Karpeev #undef __FUNCT__ 1748d709ea83SDmitry Karpeev #define __FUNCT__ "PCGASMSetDMSubdomains" 1749d709ea83SDmitry Karpeev /*@ 1750d709ea83SDmitry Karpeev PCGASMSetDMSubdomains - Indicates whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible. 1751d709ea83SDmitry Karpeev Logically Collective 1752d709ea83SDmitry Karpeev 1753d709ea83SDmitry Karpeev Input Parameter: 1754d709ea83SDmitry Karpeev + pc - the preconditioner 1755d709ea83SDmitry Karpeev - flg - boolean indicating whether to use subdomains defined by the DM 1756d709ea83SDmitry Karpeev 1757d709ea83SDmitry Karpeev Options Database Key: 1758d709ea83SDmitry Karpeev . -pc_gasm_dm_subdomains 1759d709ea83SDmitry Karpeev 1760d709ea83SDmitry Karpeev Level: intermediate 1761d709ea83SDmitry Karpeev 1762d709ea83SDmitry Karpeev Notes: 1763d709ea83SDmitry Karpeev PCGASMSetTotalSubdomains() and PCGASMSetOverlap() take precedence over PCGASMSetDMSubdomains(), 1764d709ea83SDmitry Karpeev so setting either of the first two effectively turns the latter off. 1765d709ea83SDmitry Karpeev 1766d709ea83SDmitry Karpeev .keywords: PC, ASM, DM, set, subdomains, additive Schwarz 1767d709ea83SDmitry Karpeev 1768d709ea83SDmitry Karpeev .seealso: PCGASMGetDMSubdomains(), PCGASMSetTotalSubdomains(), PCGASMSetOverlap() 1769d709ea83SDmitry Karpeev PCGASMCreateSubdomains2D() 1770d709ea83SDmitry Karpeev @*/ 1771d709ea83SDmitry Karpeev PetscErrorCode PCGASMSetDMSubdomains(PC pc,PetscBool flg) 1772d709ea83SDmitry Karpeev { 1773d709ea83SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 1774d709ea83SDmitry Karpeev PetscErrorCode ierr; 1775d709ea83SDmitry Karpeev PetscBool match; 1776d709ea83SDmitry Karpeev 1777d709ea83SDmitry Karpeev PetscFunctionBegin; 1778d709ea83SDmitry Karpeev PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1779d709ea83SDmitry Karpeev PetscValidLogicalCollectiveBool(pc,flg,2); 1780d709ea83SDmitry Karpeev if (pc->setupcalled) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Not for a setup PC."); 1781d709ea83SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match);CHKERRQ(ierr); 1782d709ea83SDmitry Karpeev if (match) { 1783d709ea83SDmitry Karpeev osm->dm_subdomains = flg; 1784d709ea83SDmitry Karpeev } 1785d709ea83SDmitry Karpeev PetscFunctionReturn(0); 1786d709ea83SDmitry Karpeev } 1787d709ea83SDmitry Karpeev 1788d709ea83SDmitry Karpeev #undef __FUNCT__ 1789d709ea83SDmitry Karpeev #define __FUNCT__ "PCGASMGetDMSubdomains" 1790d709ea83SDmitry Karpeev /*@ 1791d709ea83SDmitry Karpeev PCGASMGetDMSubdomains - Returns flag indicating whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible. 1792d709ea83SDmitry Karpeev Not Collective 1793d709ea83SDmitry Karpeev 1794d709ea83SDmitry Karpeev Input Parameter: 1795d709ea83SDmitry Karpeev . pc - the preconditioner 1796d709ea83SDmitry Karpeev 1797d709ea83SDmitry Karpeev Output Parameter: 1798d709ea83SDmitry Karpeev . flg - boolean indicating whether to use subdomains defined by the DM 1799d709ea83SDmitry Karpeev 1800d709ea83SDmitry Karpeev Level: intermediate 1801d709ea83SDmitry Karpeev 1802d709ea83SDmitry Karpeev .keywords: PC, ASM, DM, set, subdomains, additive Schwarz 1803d709ea83SDmitry Karpeev 1804d709ea83SDmitry Karpeev .seealso: PCGASMSetDMSubdomains(), PCGASMSetTotalSubdomains(), PCGASMSetOverlap() 1805d709ea83SDmitry Karpeev PCGASMCreateSubdomains2D() 1806d709ea83SDmitry Karpeev @*/ 1807d709ea83SDmitry Karpeev PetscErrorCode PCGASMGetDMSubdomains(PC pc,PetscBool* flg) 1808d709ea83SDmitry Karpeev { 1809d709ea83SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 1810d709ea83SDmitry Karpeev PetscErrorCode ierr; 1811d709ea83SDmitry Karpeev PetscBool match; 1812d709ea83SDmitry Karpeev 1813d709ea83SDmitry Karpeev PetscFunctionBegin; 1814d709ea83SDmitry Karpeev PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1815d709ea83SDmitry Karpeev PetscValidPointer(flg,2); 1816d709ea83SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match);CHKERRQ(ierr); 1817d709ea83SDmitry Karpeev if (match) { 1818d709ea83SDmitry Karpeev if (flg) *flg = osm->dm_subdomains; 1819d709ea83SDmitry Karpeev } 1820d709ea83SDmitry Karpeev PetscFunctionReturn(0); 1821d709ea83SDmitry Karpeev } 1822