1b1a0a8a3SJed Brown /* 2f746d493SDmitry Karpeev This file defines an "generalized" additive Schwarz preconditioner for any Mat implementation. 36a4f0f73SDmitry Karpeev In this version each processor may intersect multiple subdomains and any subdomain may 46a4f0f73SDmitry Karpeev intersect multiple processors. Intersections of subdomains with processors are called *local 56a4f0f73SDmitry Karpeev subdomains*. 6b1a0a8a3SJed Brown 78f3b4b4dSDmitry Karpeev N - total number of distinct global subdomains (set explicitly in PCGASMSetTotalSubdomains() or implicitly PCGASMSetSubdomains() and then calculated in PCSetUp_GASM()) 86a4f0f73SDmitry Karpeev n - actual number of local subdomains on this processor (set in PCGASMSetSubdomains() or calculated in PCGASMSetTotalSubdomains()) 98f3b4b4dSDmitry Karpeev nmax - maximum number of local subdomains per processor (calculated in PCSetUp_GASM()) 10b1a0a8a3SJed Brown */ 11af0996ceSBarry Smith #include <petsc/private/pcimpl.h> /*I "petscpc.h" I*/ 121e25c274SJed Brown #include <petscdm.h> 13b1a0a8a3SJed Brown 14b1a0a8a3SJed Brown typedef struct { 15f746d493SDmitry Karpeev PetscInt N,n,nmax; 16b1a0a8a3SJed Brown PetscInt overlap; /* overlap requested by user */ 178f3b4b4dSDmitry Karpeev PCGASMType type; /* use reduced interpolation, restriction or both */ 188f3b4b4dSDmitry Karpeev PetscBool type_set; /* if user set this value (so won't change it for symmetric problems) */ 198f3b4b4dSDmitry Karpeev PetscBool same_subdomain_solvers; /* flag indicating whether all local solvers are same */ 208f3b4b4dSDmitry Karpeev PetscBool sort_indices; /* flag to sort subdomain indices */ 218f3b4b4dSDmitry Karpeev PetscBool user_subdomains; /* whether the user set explicit subdomain index sets -- keep them on PCReset() */ 228f3b4b4dSDmitry Karpeev PetscBool dm_subdomains; /* whether DM is allowed to define subdomains */ 23ea91fabdSFande Kong PetscBool hierarchicalpartitioning; 248f3b4b4dSDmitry Karpeev IS *ois; /* index sets that define the outer (conceptually, overlapping) subdomains */ 258f3b4b4dSDmitry Karpeev IS *iis; /* index sets that define the inner (conceptually, nonoverlapping) subdomains */ 268f3b4b4dSDmitry Karpeev KSP *ksp; /* linear solvers for each subdomain */ 278f3b4b4dSDmitry Karpeev Mat *pmat; /* subdomain block matrices */ 28f746d493SDmitry Karpeev Vec gx,gy; /* Merged work vectors */ 29f746d493SDmitry Karpeev Vec *x,*y; /* Split work vectors; storage aliases pieces of storage of the above merged vectors. */ 306a4f0f73SDmitry Karpeev VecScatter gorestriction; /* merged restriction to disjoint union of outer subdomains */ 316a4f0f73SDmitry Karpeev VecScatter girestriction; /* merged restriction to disjoint union of inner subdomains */ 32ea91fabdSFande Kong VecScatter pctoouter; 33ea91fabdSFande Kong IS permutationIS; 34ea91fabdSFande Kong Mat permutationP; 35ea91fabdSFande Kong Mat pcmat; 36ea91fabdSFande Kong Vec pcx,pcy; 37f746d493SDmitry Karpeev } PC_GASM; 38b1a0a8a3SJed Brown 398f3b4b4dSDmitry Karpeev static PetscErrorCode PCGASMComputeGlobalSubdomainNumbering_Private(PC pc,PetscInt **numbering,PetscInt **permutation) 408f3b4b4dSDmitry Karpeev { 418f3b4b4dSDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 428f3b4b4dSDmitry Karpeev PetscInt i; 438f3b4b4dSDmitry Karpeev PetscErrorCode ierr; 448f3b4b4dSDmitry Karpeev 458f3b4b4dSDmitry Karpeev PetscFunctionBegin; 468f3b4b4dSDmitry Karpeev /* Determine the number of globally-distinct subdomains and compute a global numbering for them. */ 478f3b4b4dSDmitry Karpeev ierr = PetscMalloc2(osm->n,numbering,osm->n,permutation);CHKERRQ(ierr); 488f3b4b4dSDmitry Karpeev ierr = PetscObjectsListGetGlobalNumbering(PetscObjectComm((PetscObject)pc),osm->n,(PetscObject*)osm->iis,NULL,*numbering);CHKERRQ(ierr); 498f3b4b4dSDmitry Karpeev for (i = 0; i < osm->n; ++i) (*permutation)[i] = i; 508f3b4b4dSDmitry Karpeev ierr = PetscSortIntWithPermutation(osm->n,*numbering,*permutation);CHKERRQ(ierr); 518f3b4b4dSDmitry Karpeev PetscFunctionReturn(0); 528f3b4b4dSDmitry Karpeev } 538f3b4b4dSDmitry Karpeev 5406b43e7eSDmitry Karpeev static PetscErrorCode PCGASMSubdomainView_Private(PC pc, PetscInt i, PetscViewer viewer) 55af538404SDmitry Karpeev { 56af538404SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 5706b43e7eSDmitry Karpeev PetscInt j,nidx; 58af538404SDmitry Karpeev const PetscInt *idx; 5906b43e7eSDmitry Karpeev PetscViewer sviewer; 60af538404SDmitry Karpeev char *cidx; 61af538404SDmitry Karpeev PetscErrorCode ierr; 62af538404SDmitry Karpeev 63af538404SDmitry Karpeev PetscFunctionBegin; 64ce94432eSBarry 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); 654bde246dSDmitry Karpeev /* Inner subdomains. */ 664bde246dSDmitry Karpeev ierr = ISGetLocalSize(osm->iis[i], &nidx);CHKERRQ(ierr); 674bde246dSDmitry Karpeev /* 684bde246dSDmitry Karpeev No more than 15 characters per index plus a space. 694bde246dSDmitry Karpeev PetscViewerStringSPrintf requires a string of size at least 2, so use (nidx+1) instead of nidx, 704bde246dSDmitry Karpeev in case nidx == 0. That will take care of the space for the trailing '\0' as well. 714bde246dSDmitry Karpeev For nidx == 0, the whole string 16 '\0'. 724bde246dSDmitry Karpeev */ 7336a9e3b9SBarry Smith #define len 16*(nidx+1)+1 7436a9e3b9SBarry Smith ierr = PetscMalloc1(len, &cidx);CHKERRQ(ierr); 7536a9e3b9SBarry Smith ierr = PetscViewerStringOpen(PETSC_COMM_SELF, cidx, len, &sviewer);CHKERRQ(ierr); 7636a9e3b9SBarry Smith #undef len 774bde246dSDmitry Karpeev ierr = ISGetIndices(osm->iis[i], &idx);CHKERRQ(ierr); 784bde246dSDmitry Karpeev for (j = 0; j < nidx; ++j) { 794bde246dSDmitry Karpeev ierr = PetscViewerStringSPrintf(sviewer, "%D ", idx[j]);CHKERRQ(ierr); 804bde246dSDmitry Karpeev } 814bde246dSDmitry Karpeev ierr = ISRestoreIndices(osm->iis[i],&idx);CHKERRQ(ierr); 824bde246dSDmitry Karpeev ierr = PetscViewerDestroy(&sviewer);CHKERRQ(ierr); 834bde246dSDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer, "Inner subdomain:\n");CHKERRQ(ierr); 844bde246dSDmitry Karpeev ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 851575c14dSBarry Smith ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 864bde246dSDmitry Karpeev ierr = PetscViewerASCIISynchronizedPrintf(viewer, "%s", cidx);CHKERRQ(ierr); 874bde246dSDmitry Karpeev ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 881575c14dSBarry Smith ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 894bde246dSDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer, "\n");CHKERRQ(ierr); 904bde246dSDmitry Karpeev ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 914bde246dSDmitry Karpeev ierr = PetscFree(cidx);CHKERRQ(ierr); 924bde246dSDmitry Karpeev /* Outer subdomains. */ 936a4f0f73SDmitry Karpeev ierr = ISGetLocalSize(osm->ois[i], &nidx);CHKERRQ(ierr); 94eec7e3faSDmitry Karpeev /* 95eec7e3faSDmitry Karpeev No more than 15 characters per index plus a space. 96eec7e3faSDmitry Karpeev PetscViewerStringSPrintf requires a string of size at least 2, so use (nidx+1) instead of nidx, 97eec7e3faSDmitry Karpeev in case nidx == 0. That will take care of the space for the trailing '\0' as well. 98eec7e3faSDmitry Karpeev For nidx == 0, the whole string 16 '\0'. 99eec7e3faSDmitry Karpeev */ 10036a9e3b9SBarry Smith #define len 16*(nidx+1)+1 10136a9e3b9SBarry Smith ierr = PetscMalloc1(len, &cidx);CHKERRQ(ierr); 10236a9e3b9SBarry Smith ierr = PetscViewerStringOpen(PETSC_COMM_SELF, cidx, len, &sviewer);CHKERRQ(ierr); 10336a9e3b9SBarry Smith #undef len 1046a4f0f73SDmitry Karpeev ierr = ISGetIndices(osm->ois[i], &idx);CHKERRQ(ierr); 105af538404SDmitry Karpeev for (j = 0; j < nidx; ++j) { 106af538404SDmitry Karpeev ierr = PetscViewerStringSPrintf(sviewer,"%D ", idx[j]);CHKERRQ(ierr); 107af538404SDmitry Karpeev } 1086bf464f9SBarry Smith ierr = PetscViewerDestroy(&sviewer);CHKERRQ(ierr); 1096a4f0f73SDmitry Karpeev ierr = ISRestoreIndices(osm->ois[i],&idx);CHKERRQ(ierr); 1106a4f0f73SDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer, "Outer subdomain:\n");CHKERRQ(ierr); 11106b43e7eSDmitry Karpeev ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 1121575c14dSBarry Smith ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 113af538404SDmitry Karpeev ierr = PetscViewerASCIISynchronizedPrintf(viewer, "%s", cidx);CHKERRQ(ierr); 114af538404SDmitry Karpeev ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 1151575c14dSBarry Smith ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 116af538404SDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer, "\n");CHKERRQ(ierr); 11706b43e7eSDmitry Karpeev ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 118af538404SDmitry Karpeev ierr = PetscFree(cidx);CHKERRQ(ierr); 11906b43e7eSDmitry Karpeev PetscFunctionReturn(0); 12006b43e7eSDmitry Karpeev } 12106b43e7eSDmitry Karpeev 12206b43e7eSDmitry Karpeev static PetscErrorCode PCGASMPrintSubdomains(PC pc) 12306b43e7eSDmitry Karpeev { 12406b43e7eSDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 12506b43e7eSDmitry Karpeev const char *prefix; 12606b43e7eSDmitry Karpeev char fname[PETSC_MAX_PATH_LEN+1]; 1278f3b4b4dSDmitry Karpeev PetscInt l, d, count; 1288f3b4b4dSDmitry Karpeev PetscBool doprint,found; 1290298fd71SBarry Smith PetscViewer viewer, sviewer = NULL; 1308f3b4b4dSDmitry Karpeev PetscInt *numbering,*permutation;/* global numbering of locally-supported subdomains and the permutation from the local ordering */ 13106b43e7eSDmitry Karpeev PetscErrorCode ierr; 13206b43e7eSDmitry Karpeev 13306b43e7eSDmitry Karpeev PetscFunctionBegin; 13406b43e7eSDmitry Karpeev ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 1358f3b4b4dSDmitry Karpeev doprint = PETSC_FALSE; 136c5929fdfSBarry Smith ierr = PetscOptionsGetBool(NULL,prefix,"-pc_gasm_print_subdomains",&doprint,NULL);CHKERRQ(ierr); 1378f3b4b4dSDmitry Karpeev if (!doprint) PetscFunctionReturn(0); 138*589a23caSBarry Smith ierr = PetscOptionsGetString(NULL,prefix,"-pc_gasm_print_subdomains",fname,sizeof(fname),&found);CHKERRQ(ierr); 13906b43e7eSDmitry Karpeev if (!found) { ierr = PetscStrcpy(fname,"stdout");CHKERRQ(ierr); }; 140ce94432eSBarry Smith ierr = PetscViewerASCIIOpen(PetscObjectComm((PetscObject)pc),fname,&viewer);CHKERRQ(ierr); 1414bde246dSDmitry Karpeev /* 1424bde246dSDmitry Karpeev Make sure the viewer has a name. Otherwise this may cause a deadlock or other weird errors when creating a subcomm viewer: 1434bde246dSDmitry Karpeev the subcomm viewer will attempt to inherit the viewer's name, which, if not set, will be constructed collectively on the comm. 1444bde246dSDmitry Karpeev */ 1454bde246dSDmitry Karpeev ierr = PetscObjectName((PetscObject)viewer);CHKERRQ(ierr); 1464bde246dSDmitry Karpeev l = 0; 1478f3b4b4dSDmitry Karpeev ierr = PCGASMComputeGlobalSubdomainNumbering_Private(pc,&numbering,&permutation);CHKERRQ(ierr); 1488f3b4b4dSDmitry Karpeev for (count = 0; count < osm->N; ++count) { 1494bde246dSDmitry Karpeev /* Now let subdomains go one at a time in the global numbering order and print their subdomain/solver info. */ 1504bde246dSDmitry Karpeev if (l<osm->n) { 1514bde246dSDmitry Karpeev d = permutation[l]; /* d is the local number of the l-th smallest (in the global ordering) among the locally supported subdomains */ 1524bde246dSDmitry Karpeev if (numbering[d] == count) { 1533f08860eSBarry Smith ierr = PetscViewerGetSubViewer(viewer,((PetscObject)osm->ois[d])->comm, &sviewer);CHKERRQ(ierr); 1544bde246dSDmitry Karpeev ierr = PCGASMSubdomainView_Private(pc,d,sviewer);CHKERRQ(ierr); 1553f08860eSBarry Smith ierr = PetscViewerRestoreSubViewer(viewer,((PetscObject)osm->ois[d])->comm, &sviewer);CHKERRQ(ierr); 1564bde246dSDmitry Karpeev ++l; 157af538404SDmitry Karpeev } 1584bde246dSDmitry Karpeev } 159ce94432eSBarry Smith ierr = MPI_Barrier(PetscObjectComm((PetscObject)pc));CHKERRQ(ierr); 1604bde246dSDmitry Karpeev } 1618f3b4b4dSDmitry Karpeev ierr = PetscFree2(numbering,permutation);CHKERRQ(ierr); 1624bde246dSDmitry Karpeev ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 163af538404SDmitry Karpeev PetscFunctionReturn(0); 164af538404SDmitry Karpeev } 165af538404SDmitry Karpeev 166af538404SDmitry Karpeev 167f746d493SDmitry Karpeev static PetscErrorCode PCView_GASM(PC pc,PetscViewer viewer) 168b1a0a8a3SJed Brown { 169f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 170af538404SDmitry Karpeev const char *prefix; 171b1a0a8a3SJed Brown PetscErrorCode ierr; 172af538404SDmitry Karpeev PetscMPIInt rank, size; 1738f3b4b4dSDmitry Karpeev PetscInt bsz; 17406b43e7eSDmitry Karpeev PetscBool iascii,view_subdomains=PETSC_FALSE; 175b1a0a8a3SJed Brown PetscViewer sviewer; 1768f3b4b4dSDmitry Karpeev PetscInt count, l; 1776a4f0f73SDmitry Karpeev char overlap[256] = "user-defined overlap"; 1786a4f0f73SDmitry Karpeev char gsubdomains[256] = "unknown total number of subdomains"; 17906b43e7eSDmitry Karpeev char msubdomains[256] = "unknown max number of local subdomains"; 1808f3b4b4dSDmitry Karpeev PetscInt *numbering,*permutation;/* global numbering of locally-supported subdomains and the permutation from the local ordering */ 18111aeaf0aSBarry Smith 182b1a0a8a3SJed Brown PetscFunctionBegin; 183ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size);CHKERRQ(ierr); 184ce94432eSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank);CHKERRQ(ierr); 18506b43e7eSDmitry Karpeev 18606b43e7eSDmitry Karpeev if (osm->overlap >= 0) { 18706b43e7eSDmitry Karpeev ierr = PetscSNPrintf(overlap,sizeof(overlap),"requested amount of overlap = %D",osm->overlap);CHKERRQ(ierr); 18806b43e7eSDmitry Karpeev } 1898f3b4b4dSDmitry Karpeev if (osm->N != PETSC_DETERMINE) { 1908f3b4b4dSDmitry Karpeev ierr = PetscSNPrintf(gsubdomains, sizeof(gsubdomains), "total number of subdomains = %D",osm->N);CHKERRQ(ierr); 19106b43e7eSDmitry Karpeev } 1928f3b4b4dSDmitry Karpeev if (osm->nmax != PETSC_DETERMINE) { 19306b43e7eSDmitry Karpeev ierr = PetscSNPrintf(msubdomains,sizeof(msubdomains),"max number of local subdomains = %D",osm->nmax);CHKERRQ(ierr); 19406b43e7eSDmitry Karpeev } 19506b43e7eSDmitry Karpeev 196af538404SDmitry Karpeev ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 197c5929fdfSBarry Smith ierr = PetscOptionsGetBool(NULL,prefix,"-pc_gasm_view_subdomains",&view_subdomains,NULL);CHKERRQ(ierr); 19806b43e7eSDmitry Karpeev 19906b43e7eSDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 200b1a0a8a3SJed Brown if (iascii) { 20106b43e7eSDmitry Karpeev /* 20206b43e7eSDmitry Karpeev Make sure the viewer has a name. Otherwise this may cause a deadlock when creating a subcomm viewer: 20306b43e7eSDmitry Karpeev the subcomm viewer will attempt to inherit the viewer's name, which, if not set, will be constructed 20406b43e7eSDmitry Karpeev collectively on the comm. 20506b43e7eSDmitry Karpeev */ 20606b43e7eSDmitry Karpeev ierr = PetscObjectName((PetscObject)viewer);CHKERRQ(ierr); 20706b43e7eSDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer," Restriction/interpolation type: %s\n",PCGASMTypes[osm->type]);CHKERRQ(ierr); 20806b43e7eSDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer," %s\n",overlap);CHKERRQ(ierr); 20906b43e7eSDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer," %s\n",gsubdomains);CHKERRQ(ierr); 21006b43e7eSDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer," %s\n",msubdomains);CHKERRQ(ierr); 2111575c14dSBarry Smith ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 2128f3b4b4dSDmitry Karpeev ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d|%d] number of locally-supported subdomains = %D\n",rank,size,osm->n);CHKERRQ(ierr); 213af538404SDmitry Karpeev ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 2141575c14dSBarry Smith ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 2156a4f0f73SDmitry Karpeev /* Cannot take advantage of osm->same_subdomain_solvers without a global numbering of subdomains. */ 21606b43e7eSDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer," Subdomain solver info is as follows:\n");CHKERRQ(ierr); 217b1a0a8a3SJed Brown ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 218b1a0a8a3SJed Brown ierr = PetscViewerASCIIPrintf(viewer," - - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr); 21906b43e7eSDmitry Karpeev /* Make sure that everybody waits for the banner to be printed. */ 220ce94432eSBarry Smith ierr = MPI_Barrier(PetscObjectComm((PetscObject)viewer));CHKERRQ(ierr); 2214bde246dSDmitry Karpeev /* Now let subdomains go one at a time in the global numbering order and print their subdomain/solver info. */ 2228f3b4b4dSDmitry Karpeev ierr = PCGASMComputeGlobalSubdomainNumbering_Private(pc,&numbering,&permutation);CHKERRQ(ierr); 22306b43e7eSDmitry Karpeev l = 0; 2248f3b4b4dSDmitry Karpeev for (count = 0; count < osm->N; ++count) { 22506b43e7eSDmitry Karpeev PetscMPIInt srank, ssize; 22606b43e7eSDmitry Karpeev if (l<osm->n) { 22706b43e7eSDmitry Karpeev PetscInt d = permutation[l]; /* d is the local number of the l-th smallest (in the global ordering) among the locally supported subdomains */ 22806b43e7eSDmitry Karpeev if (numbering[d] == count) { 2296a4f0f73SDmitry Karpeev ierr = MPI_Comm_size(((PetscObject)osm->ois[d])->comm, &ssize);CHKERRQ(ierr); 2306a4f0f73SDmitry Karpeev ierr = MPI_Comm_rank(((PetscObject)osm->ois[d])->comm, &srank);CHKERRQ(ierr); 2313f08860eSBarry Smith ierr = PetscViewerGetSubViewer(viewer,((PetscObject)osm->ois[d])->comm, &sviewer);CHKERRQ(ierr); 2326a4f0f73SDmitry Karpeev ierr = ISGetLocalSize(osm->ois[d],&bsz);CHKERRQ(ierr); 2338f3b4b4dSDmitry Karpeev ierr = PetscViewerASCIISynchronizedPrintf(sviewer," [%d|%d] (subcomm [%d|%d]) local subdomain number %D, local size = %D\n",rank,size,srank,ssize,d,bsz);CHKERRQ(ierr); 2346a4f0f73SDmitry Karpeev ierr = PetscViewerFlush(sviewer);CHKERRQ(ierr); 23506b43e7eSDmitry Karpeev if (view_subdomains) { 23606b43e7eSDmitry Karpeev ierr = PCGASMSubdomainView_Private(pc,d,sviewer);CHKERRQ(ierr); 23706b43e7eSDmitry Karpeev } 2386a4f0f73SDmitry Karpeev if (!pc->setupcalled) { 2391575c14dSBarry Smith ierr = PetscViewerASCIIPrintf(sviewer, " Solver not set up yet: PCSetUp() not yet called\n");CHKERRQ(ierr); 2402fa5cd67SKarl Rupp } else { 24106b43e7eSDmitry Karpeev ierr = KSPView(osm->ksp[d],sviewer);CHKERRQ(ierr); 2426a4f0f73SDmitry Karpeev } 24306b43e7eSDmitry Karpeev ierr = PetscViewerASCIIPrintf(sviewer," - - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr); 2447b23a99aSBarry Smith ierr = PetscViewerFlush(sviewer);CHKERRQ(ierr); 2453f08860eSBarry Smith ierr = PetscViewerRestoreSubViewer(viewer,((PetscObject)osm->ois[d])->comm, &sviewer);CHKERRQ(ierr); 24606b43e7eSDmitry Karpeev ++l; 247b1a0a8a3SJed Brown } 24806b43e7eSDmitry Karpeev } 249ce94432eSBarry Smith ierr = MPI_Barrier(PetscObjectComm((PetscObject)pc));CHKERRQ(ierr); 250b1a0a8a3SJed Brown } 2518f3b4b4dSDmitry Karpeev ierr = PetscFree2(numbering,permutation);CHKERRQ(ierr); 252b1a0a8a3SJed Brown ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 253b1a0a8a3SJed Brown } 254c5e4d11fSDmitry Karpeev ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 2551575c14dSBarry Smith ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 256b1a0a8a3SJed Brown PetscFunctionReturn(0); 257b1a0a8a3SJed Brown } 258b1a0a8a3SJed Brown 2598f3b4b4dSDmitry Karpeev PETSC_INTERN PetscErrorCode PCGASMCreateLocalSubdomains(Mat A, PetscInt nloc, IS *iis[]); 260af538404SDmitry Karpeev 261ea91fabdSFande Kong 262ea91fabdSFande Kong 263ea91fabdSFande Kong PetscErrorCode PCGASMSetHierarchicalPartitioning(PC pc) 264ea91fabdSFande Kong { 265ea91fabdSFande Kong PC_GASM *osm = (PC_GASM*)pc->data; 266ea91fabdSFande Kong MatPartitioning part; 267ea91fabdSFande Kong MPI_Comm comm; 268ea91fabdSFande Kong PetscMPIInt size; 269ea91fabdSFande Kong PetscInt nlocalsubdomains,fromrows_localsize; 270ea91fabdSFande Kong IS partitioning,fromrows,isn; 271ea91fabdSFande Kong Vec outervec; 272ea91fabdSFande Kong PetscErrorCode ierr; 273ea91fabdSFande Kong 274ea91fabdSFande Kong PetscFunctionBegin; 275ea91fabdSFande Kong ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 276ea91fabdSFande Kong ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 277ea91fabdSFande Kong /* we do not need a hierarchical partitioning when 278ea91fabdSFande Kong * the total number of subdomains is consistent with 279ea91fabdSFande Kong * the number of MPI tasks. 280ea91fabdSFande Kong * For the following cases, we do not need to use HP 281ea91fabdSFande Kong * */ 282670417b2SFande Kong if(osm->N==PETSC_DETERMINE || osm->N>=size || osm->N==1) PetscFunctionReturn(0); 283670417b2SFande Kong if(size%osm->N != 0) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_ARG_INCOMP,"have to specify the total number of subdomains %D to be a factor of the number of processors %d \n",osm->N,size); 284ea91fabdSFande Kong nlocalsubdomains = size/osm->N; 285ea91fabdSFande Kong osm->n = 1; 286ea91fabdSFande Kong ierr = MatPartitioningCreate(comm,&part);CHKERRQ(ierr); 287ea91fabdSFande Kong ierr = MatPartitioningSetAdjacency(part,pc->pmat);CHKERRQ(ierr); 288ea91fabdSFande Kong ierr = MatPartitioningSetType(part,MATPARTITIONINGHIERARCH);CHKERRQ(ierr); 289ea91fabdSFande Kong ierr = MatPartitioningHierarchicalSetNcoarseparts(part,osm->N);CHKERRQ(ierr); 290ea91fabdSFande Kong ierr = MatPartitioningHierarchicalSetNfineparts(part,nlocalsubdomains);CHKERRQ(ierr); 291ea91fabdSFande Kong ierr = MatPartitioningSetFromOptions(part);CHKERRQ(ierr); 292ea91fabdSFande Kong /* get new processor owner number of each vertex */ 293ea91fabdSFande Kong ierr = MatPartitioningApply(part,&partitioning);CHKERRQ(ierr); 294ea91fabdSFande Kong ierr = ISBuildTwoSided(partitioning,NULL,&fromrows);CHKERRQ(ierr); 295ea91fabdSFande Kong ierr = ISPartitioningToNumbering(partitioning,&isn);CHKERRQ(ierr); 296ea91fabdSFande Kong ierr = ISDestroy(&isn);CHKERRQ(ierr); 297ea91fabdSFande Kong ierr = ISGetLocalSize(fromrows,&fromrows_localsize);CHKERRQ(ierr); 298ea91fabdSFande Kong ierr = MatPartitioningDestroy(&part);CHKERRQ(ierr); 299ea91fabdSFande Kong ierr = MatCreateVecs(pc->pmat,&outervec,NULL);CHKERRQ(ierr); 300ea91fabdSFande Kong ierr = VecCreateMPI(comm,fromrows_localsize,PETSC_DETERMINE,&(osm->pcx));CHKERRQ(ierr); 301ea91fabdSFande Kong ierr = VecDuplicate(osm->pcx,&(osm->pcy));CHKERRQ(ierr); 3029448b7f1SJunchao Zhang ierr = VecScatterCreate(osm->pcx,NULL,outervec,fromrows,&(osm->pctoouter));CHKERRQ(ierr); 3037dae84e0SHong Zhang ierr = MatCreateSubMatrix(pc->pmat,fromrows,fromrows,MAT_INITIAL_MATRIX,&(osm->permutationP));CHKERRQ(ierr); 304ea91fabdSFande Kong ierr = PetscObjectReference((PetscObject)fromrows);CHKERRQ(ierr); 305ea91fabdSFande Kong osm->permutationIS = fromrows; 306ea91fabdSFande Kong osm->pcmat = pc->pmat; 307ea91fabdSFande Kong ierr = PetscObjectReference((PetscObject)osm->permutationP);CHKERRQ(ierr); 308ea91fabdSFande Kong pc->pmat = osm->permutationP; 309ea91fabdSFande Kong ierr = VecDestroy(&outervec);CHKERRQ(ierr); 310ea91fabdSFande Kong ierr = ISDestroy(&fromrows);CHKERRQ(ierr); 311ea91fabdSFande Kong ierr = ISDestroy(&partitioning);CHKERRQ(ierr); 312ea91fabdSFande Kong osm->n = PETSC_DETERMINE; 313ea91fabdSFande Kong PetscFunctionReturn(0); 314ea91fabdSFande Kong } 315ea91fabdSFande Kong 316ea91fabdSFande Kong 317ea91fabdSFande Kong 318f746d493SDmitry Karpeev static PetscErrorCode PCSetUp_GASM(PC pc) 319b1a0a8a3SJed Brown { 320f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 321b1a0a8a3SJed Brown PetscErrorCode ierr; 322ea91fabdSFande Kong PetscInt i,nInnerIndices,nTotalInnerIndices; 323c10bc1a1SDmitry Karpeev PetscMPIInt rank, size; 324b1a0a8a3SJed Brown MatReuse scall = MAT_REUSE_MATRIX; 325b1a0a8a3SJed Brown KSP ksp; 326b1a0a8a3SJed Brown PC subpc; 327b1a0a8a3SJed Brown const char *prefix,*pprefix; 328f746d493SDmitry Karpeev Vec x,y; 3296a4f0f73SDmitry Karpeev PetscInt oni; /* Number of indices in the i-th local outer subdomain. */ 3306a4f0f73SDmitry Karpeev const PetscInt *oidxi; /* Indices from the i-th subdomain local outer subdomain. */ 3316a4f0f73SDmitry Karpeev PetscInt on; /* Number of indices in the disjoint union of local outer subdomains. */ 3326a4f0f73SDmitry Karpeev PetscInt *oidx; /* Indices in the disjoint union of local outer subdomains. */ 3336a4f0f73SDmitry Karpeev IS gois; /* Disjoint union the global indices of outer subdomains. */ 3346a4f0f73SDmitry Karpeev IS goid; /* Identity IS of the size of the disjoint union of outer subdomains. */ 335f746d493SDmitry Karpeev PetscScalar *gxarray, *gyarray; 336930d09c1SFande Kong PetscInt gostart; /* Start of locally-owned indices in the vectors -- osm->gx,osm->gy -- over the disjoint union of outer subdomains. */ 3378f3b4b4dSDmitry Karpeev PetscInt num_subdomains = 0; 3380298fd71SBarry Smith DM *subdomain_dm = NULL; 3398f3b4b4dSDmitry Karpeev char **subdomain_names = NULL; 340f771a274SFande Kong PetscInt *numbering; 3418f3b4b4dSDmitry Karpeev 342b1a0a8a3SJed Brown 343b1a0a8a3SJed Brown PetscFunctionBegin; 344ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr); 345ce94432eSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);CHKERRQ(ierr); 346b1a0a8a3SJed Brown if (!pc->setupcalled) { 347ea91fabdSFande Kong /* use a hierarchical partitioning */ 348ea91fabdSFande Kong if(osm->hierarchicalpartitioning){ 349ea91fabdSFande Kong ierr = PCGASMSetHierarchicalPartitioning(pc);CHKERRQ(ierr); 350ea91fabdSFande Kong } 3518f3b4b4dSDmitry Karpeev if (osm->n == PETSC_DETERMINE) { 3528f3b4b4dSDmitry Karpeev if (osm->N != PETSC_DETERMINE) { 3538f3b4b4dSDmitry Karpeev /* No local subdomains given, but the desired number of total subdomains is known, so construct them accordingly. */ 3548f3b4b4dSDmitry Karpeev ierr = PCGASMCreateSubdomains(pc->pmat,osm->N,&osm->n,&osm->iis);CHKERRQ(ierr); 3558f3b4b4dSDmitry Karpeev } else if (osm->dm_subdomains && pc->dm) { 3568f3b4b4dSDmitry Karpeev /* try pc->dm next, if allowed */ 3578f3b4b4dSDmitry Karpeev PetscInt d; 358a35b7b57SDmitry Karpeev IS *inner_subdomain_is, *outer_subdomain_is; 359a35b7b57SDmitry Karpeev ierr = DMCreateDomainDecomposition(pc->dm, &num_subdomains, &subdomain_names, &inner_subdomain_is, &outer_subdomain_is, &subdomain_dm);CHKERRQ(ierr); 360a35b7b57SDmitry Karpeev if (num_subdomains) { 361a35b7b57SDmitry Karpeev ierr = PCGASMSetSubdomains(pc, num_subdomains, inner_subdomain_is, outer_subdomain_is);CHKERRQ(ierr); 36269ca1f37SDmitry Karpeev } 363a35b7b57SDmitry Karpeev for (d = 0; d < num_subdomains; ++d) { 364a35b7b57SDmitry Karpeev if (inner_subdomain_is) {ierr = ISDestroy(&inner_subdomain_is[d]);CHKERRQ(ierr);} 365a35b7b57SDmitry Karpeev if (outer_subdomain_is) {ierr = ISDestroy(&outer_subdomain_is[d]);CHKERRQ(ierr);} 366fdc48646SDmitry Karpeev } 367a35b7b57SDmitry Karpeev ierr = PetscFree(inner_subdomain_is);CHKERRQ(ierr); 368a35b7b57SDmitry Karpeev ierr = PetscFree(outer_subdomain_is);CHKERRQ(ierr); 3698f3b4b4dSDmitry Karpeev } else { 3708f3b4b4dSDmitry Karpeev /* still no subdomains; use one per processor */ 37144fe18e5SDmitry Karpeev osm->nmax = osm->n = 1; 372ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr); 373f746d493SDmitry Karpeev osm->N = size; 3748f3b4b4dSDmitry Karpeev ierr = PCGASMCreateLocalSubdomains(pc->pmat,osm->n,&osm->iis);CHKERRQ(ierr); 375fdc48646SDmitry Karpeev } 37606b43e7eSDmitry Karpeev } 377a35b7b57SDmitry Karpeev if (!osm->iis) { 378a35b7b57SDmitry Karpeev /* 3798f3b4b4dSDmitry Karpeev osm->n was set in PCGASMSetSubdomains(), but the actual subdomains have not been supplied. 3808f3b4b4dSDmitry Karpeev We create the requisite number of local inner subdomains and then expand them into 3818f3b4b4dSDmitry Karpeev out subdomains, if necessary. 382a35b7b57SDmitry Karpeev */ 3838f3b4b4dSDmitry Karpeev ierr = PCGASMCreateLocalSubdomains(pc->pmat,osm->n,&osm->iis);CHKERRQ(ierr); 384a35b7b57SDmitry Karpeev } 3858f3b4b4dSDmitry Karpeev if (!osm->ois) { 3868f3b4b4dSDmitry Karpeev /* 3878f3b4b4dSDmitry Karpeev Initially make outer subdomains the same as inner subdomains. If nonzero additional overlap 3888f3b4b4dSDmitry Karpeev has been requested, copy the inner subdomains over so they can be modified. 3898f3b4b4dSDmitry Karpeev */ 3908f3b4b4dSDmitry Karpeev ierr = PetscMalloc1(osm->n,&osm->ois);CHKERRQ(ierr); 3918f3b4b4dSDmitry Karpeev for (i=0; i<osm->n; ++i) { 392ea91fabdSFande Kong if (osm->overlap > 0 && osm->N>1) { /* With positive overlap, osm->iis[i] will be modified */ 3938f3b4b4dSDmitry Karpeev ierr = ISDuplicate(osm->iis[i],(osm->ois)+i);CHKERRQ(ierr); 3948f3b4b4dSDmitry Karpeev ierr = ISCopy(osm->iis[i],osm->ois[i]);CHKERRQ(ierr); 3958f3b4b4dSDmitry Karpeev } else { 3968f3b4b4dSDmitry Karpeev ierr = PetscObjectReference((PetscObject)((osm->iis)[i]));CHKERRQ(ierr); 3978f3b4b4dSDmitry Karpeev osm->ois[i] = osm->iis[i]; 3988f3b4b4dSDmitry Karpeev } 3998f3b4b4dSDmitry Karpeev } 400ea91fabdSFande Kong if (osm->overlap>0 && osm->N>1) { 4018f3b4b4dSDmitry Karpeev /* Extend the "overlapping" regions by a number of steps */ 402d21f2a47SFande Kong ierr = MatIncreaseOverlapSplit(pc->pmat,osm->n,osm->ois,osm->overlap);CHKERRQ(ierr); 4038f3b4b4dSDmitry Karpeev } 404b1a0a8a3SJed Brown } 4056a4f0f73SDmitry Karpeev 4068f3b4b4dSDmitry Karpeev /* Now the subdomains are defined. Determine their global and max local numbers, if necessary. */ 4078f3b4b4dSDmitry Karpeev if (osm->nmax == PETSC_DETERMINE) { 4088f3b4b4dSDmitry Karpeev PetscMPIInt inwork,outwork; 4098f3b4b4dSDmitry Karpeev /* determine global number of subdomains and the max number of local subdomains */ 4108f3b4b4dSDmitry Karpeev inwork = osm->n; 411b2566f29SBarry Smith ierr = MPIU_Allreduce(&inwork,&outwork,1,MPI_INT,MPI_MAX,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr); 4128f3b4b4dSDmitry Karpeev osm->nmax = outwork; 4138f3b4b4dSDmitry Karpeev } 4148f3b4b4dSDmitry Karpeev if (osm->N == PETSC_DETERMINE) { 4158f3b4b4dSDmitry Karpeev /* Determine the number of globally-distinct subdomains and compute a global numbering for them. */ 4168f3b4b4dSDmitry Karpeev ierr = PetscObjectsListGetGlobalNumbering(PetscObjectComm((PetscObject)pc),osm->n,(PetscObject*)osm->ois,&osm->N,NULL);CHKERRQ(ierr); 4178f3b4b4dSDmitry Karpeev } 4188f3b4b4dSDmitry Karpeev 419b1a0a8a3SJed Brown 420b1a0a8a3SJed Brown if (osm->sort_indices) { 421f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 4226a4f0f73SDmitry Karpeev ierr = ISSort(osm->ois[i]);CHKERRQ(ierr); 4236a4f0f73SDmitry Karpeev ierr = ISSort(osm->iis[i]);CHKERRQ(ierr); 424b1a0a8a3SJed Brown } 425b1a0a8a3SJed Brown } 4268f3b4b4dSDmitry Karpeev ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 4278f3b4b4dSDmitry Karpeev ierr = PCGASMPrintSubdomains(pc);CHKERRQ(ierr); 4288f3b4b4dSDmitry Karpeev 4296a4f0f73SDmitry Karpeev /* 4306a4f0f73SDmitry Karpeev Merge the ISs, create merged vectors and restrictions. 4316a4f0f73SDmitry Karpeev */ 4326a4f0f73SDmitry Karpeev /* Merge outer subdomain ISs and construct a restriction onto the disjoint union of local outer subdomains. */ 4336a4f0f73SDmitry Karpeev on = 0; 434f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 4356a4f0f73SDmitry Karpeev ierr = ISGetLocalSize(osm->ois[i],&oni);CHKERRQ(ierr); 4366a4f0f73SDmitry Karpeev on += oni; 437f746d493SDmitry Karpeev } 438785e854fSJed Brown ierr = PetscMalloc1(on, &oidx);CHKERRQ(ierr); 4396a4f0f73SDmitry Karpeev on = 0; 440930d09c1SFande Kong /* Merge local indices together */ 441f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 4426a4f0f73SDmitry Karpeev ierr = ISGetLocalSize(osm->ois[i],&oni);CHKERRQ(ierr); 4436a4f0f73SDmitry Karpeev ierr = ISGetIndices(osm->ois[i],&oidxi);CHKERRQ(ierr); 444580bdb30SBarry Smith ierr = PetscArraycpy(oidx+on,oidxi,oni);CHKERRQ(ierr); 4456a4f0f73SDmitry Karpeev ierr = ISRestoreIndices(osm->ois[i],&oidxi);CHKERRQ(ierr); 4466a4f0f73SDmitry Karpeev on += oni; 447f746d493SDmitry Karpeev } 4486a4f0f73SDmitry Karpeev ierr = ISCreateGeneral(((PetscObject)(pc))->comm,on,oidx,PETSC_OWN_POINTER,&gois);CHKERRQ(ierr); 449ea91fabdSFande Kong nTotalInnerIndices = 0; 450ea91fabdSFande Kong for(i=0; i<osm->n; i++){ 451ea91fabdSFande Kong ierr = ISGetLocalSize(osm->iis[i],&nInnerIndices);CHKERRQ(ierr); 452ea91fabdSFande Kong nTotalInnerIndices += nInnerIndices; 453ea91fabdSFande Kong } 454ea91fabdSFande Kong ierr = VecCreateMPI(((PetscObject)(pc))->comm,nTotalInnerIndices,PETSC_DETERMINE,&x);CHKERRQ(ierr); 455ea91fabdSFande Kong ierr = VecDuplicate(x,&y);CHKERRQ(ierr); 456ea91fabdSFande Kong 457ce94432eSBarry Smith ierr = VecCreateMPI(PetscObjectComm((PetscObject)pc),on,PETSC_DECIDE,&osm->gx);CHKERRQ(ierr); 458f746d493SDmitry Karpeev ierr = VecDuplicate(osm->gx,&osm->gy);CHKERRQ(ierr); 459930d09c1SFande Kong ierr = VecGetOwnershipRange(osm->gx, &gostart, NULL);CHKERRQ(ierr); 460930d09c1SFande Kong ierr = ISCreateStride(PetscObjectComm((PetscObject)pc),on,gostart,1, &goid);CHKERRQ(ierr); 461930d09c1SFande Kong /* gois might indices not on local */ 4629448b7f1SJunchao Zhang ierr = VecScatterCreate(x,gois,osm->gx,goid, &(osm->gorestriction));CHKERRQ(ierr); 463f771a274SFande Kong ierr = PetscMalloc1(osm->n,&numbering);CHKERRQ(ierr); 464f771a274SFande Kong ierr = PetscObjectsListGetGlobalNumbering(PetscObjectComm((PetscObject)pc),osm->n,(PetscObject*)osm->ois,NULL,numbering);CHKERRQ(ierr); 4656a4f0f73SDmitry Karpeev ierr = VecDestroy(&x);CHKERRQ(ierr); 4667105d80bSDmitry Karpeev ierr = ISDestroy(&gois);CHKERRQ(ierr); 4672fa5cd67SKarl Rupp 4686a4f0f73SDmitry Karpeev /* Merge inner subdomain ISs and construct a restriction onto the disjoint union of local inner subdomains. */ 4692fa5cd67SKarl Rupp { 4702fa5cd67SKarl Rupp PetscInt ini; /* Number of indices the i-th a local inner subdomain. */ 4718966356dSPierre Jolivet PetscInt in; /* Number of indices in the disjoint union of local inner subdomains. */ 4726a4f0f73SDmitry Karpeev PetscInt *iidx; /* Global indices in the merged local inner subdomain. */ 4736a4f0f73SDmitry Karpeev PetscInt *ioidx; /* Global indices of the disjoint union of inner subdomains within the disjoint union of outer subdomains. */ 4746a4f0f73SDmitry Karpeev IS giis; /* IS for the disjoint union of inner subdomains. */ 4756a4f0f73SDmitry Karpeev IS giois; /* IS for the disjoint union of inner subdomains within the disjoint union of outer subdomains. */ 476f771a274SFande Kong PetscScalar *array; 477f771a274SFande Kong const PetscInt *indices; 478f771a274SFande Kong PetscInt k; 4796a4f0f73SDmitry Karpeev on = 0; 480f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 4816a4f0f73SDmitry Karpeev ierr = ISGetLocalSize(osm->ois[i],&oni);CHKERRQ(ierr); 4826a4f0f73SDmitry Karpeev on += oni; 483f746d493SDmitry Karpeev } 484f771a274SFande Kong ierr = PetscMalloc1(on, &iidx);CHKERRQ(ierr); 485f771a274SFande Kong ierr = PetscMalloc1(on, &ioidx);CHKERRQ(ierr); 486f771a274SFande Kong ierr = VecGetArray(y,&array);CHKERRQ(ierr); 487e12b4729SFande Kong /* set communicator id to determine where overlap is */ 488f771a274SFande Kong in = 0; 489f771a274SFande Kong for (i=0; i<osm->n; i++) { 490f771a274SFande Kong ierr = ISGetLocalSize(osm->iis[i],&ini);CHKERRQ(ierr); 491f771a274SFande Kong for (k = 0; k < ini; ++k){ 492f771a274SFande Kong array[in+k] = numbering[i]; 493f771a274SFande Kong } 494f771a274SFande Kong in += ini; 495f771a274SFande Kong } 496f771a274SFande Kong ierr = VecRestoreArray(y,&array);CHKERRQ(ierr); 497f771a274SFande Kong ierr = VecScatterBegin(osm->gorestriction,y,osm->gy,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 498f771a274SFande Kong ierr = VecScatterEnd(osm->gorestriction,y,osm->gy,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 499f771a274SFande Kong ierr = VecGetOwnershipRange(osm->gy,&gostart, NULL);CHKERRQ(ierr); 500f771a274SFande Kong ierr = VecGetArray(osm->gy,&array);CHKERRQ(ierr); 501f771a274SFande Kong on = 0; 502f771a274SFande Kong in = 0; 503f771a274SFande Kong for (i=0; i<osm->n; i++) { 504f771a274SFande Kong ierr = ISGetLocalSize(osm->ois[i],&oni);CHKERRQ(ierr); 505f771a274SFande Kong ierr = ISGetIndices(osm->ois[i],&indices);CHKERRQ(ierr); 506f771a274SFande Kong for (k=0; k<oni; k++) { 507e12b4729SFande Kong /* skip overlapping indices to get inner domain */ 50843081b6cSDmitry Karpeev if(PetscRealPart(array[on+k]) != numbering[i]) continue; 509f771a274SFande Kong iidx[in] = indices[k]; 510f771a274SFande Kong ioidx[in++] = gostart+on+k; 511f771a274SFande Kong } 512f771a274SFande Kong ierr = ISRestoreIndices(osm->ois[i], &indices);CHKERRQ(ierr); 513f771a274SFande Kong on += oni; 514f771a274SFande Kong } 515f771a274SFande Kong ierr = VecRestoreArray(osm->gy,&array);CHKERRQ(ierr); 516ce94432eSBarry Smith ierr = ISCreateGeneral(PetscObjectComm((PetscObject)pc),in,iidx,PETSC_OWN_POINTER,&giis);CHKERRQ(ierr); 517ce94432eSBarry Smith ierr = ISCreateGeneral(PetscObjectComm((PetscObject)pc),in,ioidx,PETSC_OWN_POINTER,&giois);CHKERRQ(ierr); 5189448b7f1SJunchao Zhang ierr = VecScatterCreate(y,giis,osm->gy,giois,&osm->girestriction);CHKERRQ(ierr); 5196a4f0f73SDmitry Karpeev ierr = VecDestroy(&y);CHKERRQ(ierr); 5206a4f0f73SDmitry Karpeev ierr = ISDestroy(&giis);CHKERRQ(ierr); 5216a4f0f73SDmitry Karpeev ierr = ISDestroy(&giois);CHKERRQ(ierr); 522b1a0a8a3SJed Brown } 5236a4f0f73SDmitry Karpeev ierr = ISDestroy(&goid);CHKERRQ(ierr); 524f771a274SFande Kong ierr = PetscFree(numbering);CHKERRQ(ierr); 5252fa5cd67SKarl Rupp 526f746d493SDmitry Karpeev /* Create the subdomain work vectors. */ 527785e854fSJed Brown ierr = PetscMalloc1(osm->n,&osm->x);CHKERRQ(ierr); 528785e854fSJed Brown ierr = PetscMalloc1(osm->n,&osm->y);CHKERRQ(ierr); 529f746d493SDmitry Karpeev ierr = VecGetArray(osm->gx, &gxarray);CHKERRQ(ierr); 530f746d493SDmitry Karpeev ierr = VecGetArray(osm->gy, &gyarray);CHKERRQ(ierr); 5316a4f0f73SDmitry Karpeev for (i=0, on=0; i<osm->n; ++i, on += oni) { 5326a4f0f73SDmitry Karpeev PetscInt oNi; 5336a4f0f73SDmitry Karpeev ierr = ISGetLocalSize(osm->ois[i],&oni);CHKERRQ(ierr); 534930d09c1SFande Kong /* on a sub communicator */ 5356a4f0f73SDmitry Karpeev ierr = ISGetSize(osm->ois[i],&oNi);CHKERRQ(ierr); 5366a4f0f73SDmitry Karpeev ierr = VecCreateMPIWithArray(((PetscObject)(osm->ois[i]))->comm,1,oni,oNi,gxarray+on,&osm->x[i]);CHKERRQ(ierr); 5376a4f0f73SDmitry Karpeev ierr = VecCreateMPIWithArray(((PetscObject)(osm->ois[i]))->comm,1,oni,oNi,gyarray+on,&osm->y[i]);CHKERRQ(ierr); 538b1a0a8a3SJed Brown } 539f746d493SDmitry Karpeev ierr = VecRestoreArray(osm->gx, &gxarray);CHKERRQ(ierr); 540f746d493SDmitry Karpeev ierr = VecRestoreArray(osm->gy, &gyarray);CHKERRQ(ierr); 5418f3b4b4dSDmitry Karpeev /* Create the subdomain solvers */ 542785e854fSJed Brown ierr = PetscMalloc1(osm->n,&osm->ksp);CHKERRQ(ierr); 5438f3b4b4dSDmitry Karpeev for (i=0; i<osm->n; i++) { 5448f3b4b4dSDmitry Karpeev char subprefix[PETSC_MAX_PATH_LEN+1]; 5456a4f0f73SDmitry Karpeev ierr = KSPCreate(((PetscObject)(osm->ois[i]))->comm,&ksp);CHKERRQ(ierr); 54666b14c46SBarry Smith ierr = KSPSetErrorIfNotConverged(ksp,pc->erroriffailure);CHKERRQ(ierr); 5473bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ksp);CHKERRQ(ierr); 548b1a0a8a3SJed Brown ierr = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,1);CHKERRQ(ierr); 549b1a0a8a3SJed Brown ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr); 5508f3b4b4dSDmitry Karpeev ierr = KSPGetPC(ksp,&subpc);CHKERRQ(ierr); /* Why do we need this here? */ 5518f3b4b4dSDmitry Karpeev if (subdomain_dm) { 5528f3b4b4dSDmitry Karpeev ierr = KSPSetDM(ksp,subdomain_dm[i]);CHKERRQ(ierr); 5538f3b4b4dSDmitry Karpeev ierr = DMDestroy(subdomain_dm+i);CHKERRQ(ierr); 5548f3b4b4dSDmitry Karpeev } 555b1a0a8a3SJed Brown ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 556b1a0a8a3SJed Brown ierr = KSPSetOptionsPrefix(ksp,prefix);CHKERRQ(ierr); 5578f3b4b4dSDmitry Karpeev if (subdomain_names && subdomain_names[i]) { 5588f3b4b4dSDmitry Karpeev ierr = PetscSNPrintf(subprefix,PETSC_MAX_PATH_LEN,"sub_%s_",subdomain_names[i]);CHKERRQ(ierr); 5598f3b4b4dSDmitry Karpeev ierr = KSPAppendOptionsPrefix(ksp,subprefix);CHKERRQ(ierr); 5608f3b4b4dSDmitry Karpeev ierr = PetscFree(subdomain_names[i]);CHKERRQ(ierr); 5618f3b4b4dSDmitry Karpeev } 562b1a0a8a3SJed Brown ierr = KSPAppendOptionsPrefix(ksp,"sub_");CHKERRQ(ierr); 563b1a0a8a3SJed Brown osm->ksp[i] = ksp; 564b1a0a8a3SJed Brown } 5658f3b4b4dSDmitry Karpeev ierr = PetscFree(subdomain_dm);CHKERRQ(ierr); 5668f3b4b4dSDmitry Karpeev ierr = PetscFree(subdomain_names);CHKERRQ(ierr); 567b1a0a8a3SJed Brown scall = MAT_INITIAL_MATRIX; 568b1a0a8a3SJed Brown 5698f3b4b4dSDmitry Karpeev } else { /* if (pc->setupcalled) */ 570b1a0a8a3SJed Brown /* 5718f3b4b4dSDmitry Karpeev Destroy the submatrices from the previous iteration 572b1a0a8a3SJed Brown */ 573b1a0a8a3SJed Brown if (pc->flag == DIFFERENT_NONZERO_PATTERN) { 574f746d493SDmitry Karpeev ierr = MatDestroyMatrices(osm->n,&osm->pmat);CHKERRQ(ierr); 575b1a0a8a3SJed Brown scall = MAT_INITIAL_MATRIX; 576b1a0a8a3SJed Brown } 577ea91fabdSFande Kong if(osm->permutationIS){ 5787dae84e0SHong Zhang ierr = MatCreateSubMatrix(pc->pmat,osm->permutationIS,osm->permutationIS,scall,&osm->permutationP);CHKERRQ(ierr); 579ea91fabdSFande Kong ierr = PetscObjectReference((PetscObject)osm->permutationP);CHKERRQ(ierr); 580ea91fabdSFande Kong osm->pcmat = pc->pmat; 581ea91fabdSFande Kong pc->pmat = osm->permutationP; 582b1a0a8a3SJed Brown } 583b1a0a8a3SJed Brown 584ea91fabdSFande Kong } 585ea91fabdSFande Kong 586ea91fabdSFande Kong 587b1a0a8a3SJed Brown /* 588f746d493SDmitry Karpeev Extract out the submatrices. 589b1a0a8a3SJed Brown */ 59081a5c4aaSDmitry Karpeev if (size > 1) { 5917dae84e0SHong Zhang ierr = MatCreateSubMatricesMPI(pc->pmat,osm->n,osm->ois,osm->ois,scall,&osm->pmat);CHKERRQ(ierr); 5922fa5cd67SKarl Rupp } else { 5937dae84e0SHong Zhang ierr = MatCreateSubMatrices(pc->pmat,osm->n,osm->ois,osm->ois,scall,&osm->pmat);CHKERRQ(ierr); 59481a5c4aaSDmitry Karpeev } 595b1a0a8a3SJed Brown if (scall == MAT_INITIAL_MATRIX) { 596b1a0a8a3SJed Brown ierr = PetscObjectGetOptionsPrefix((PetscObject)pc->pmat,&pprefix);CHKERRQ(ierr); 597f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 5983bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)osm->pmat[i]);CHKERRQ(ierr); 599b1a0a8a3SJed Brown ierr = PetscObjectSetOptionsPrefix((PetscObject)osm->pmat[i],pprefix);CHKERRQ(ierr); 600b1a0a8a3SJed Brown } 601b1a0a8a3SJed Brown } 602b1a0a8a3SJed Brown 603b1a0a8a3SJed Brown /* Return control to the user so that the submatrices can be modified (e.g., to apply 604b1a0a8a3SJed Brown different boundary conditions for the submatrices than for the global problem) */ 6056a4f0f73SDmitry Karpeev ierr = PCModifySubMatrices(pc,osm->n,osm->ois,osm->ois,osm->pmat,pc->modifysubmatricesP);CHKERRQ(ierr); 606b1a0a8a3SJed Brown 607b1a0a8a3SJed Brown /* 6086a4f0f73SDmitry Karpeev Loop over submatrices putting them into local ksps 609b1a0a8a3SJed Brown */ 610f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 61123ee1639SBarry Smith ierr = KSPSetOperators(osm->ksp[i],osm->pmat[i],osm->pmat[i]);CHKERRQ(ierr); 612b1a0a8a3SJed Brown if (!pc->setupcalled) { 613b1a0a8a3SJed Brown ierr = KSPSetFromOptions(osm->ksp[i]);CHKERRQ(ierr); 614b1a0a8a3SJed Brown } 615b1a0a8a3SJed Brown } 616ea91fabdSFande Kong if(osm->pcmat){ 617ea91fabdSFande Kong ierr = MatDestroy(&pc->pmat);CHKERRQ(ierr); 618ea91fabdSFande Kong pc->pmat = osm->pcmat; 6190a545947SLisandro Dalcin osm->pcmat = NULL; 620ea91fabdSFande Kong } 621b1a0a8a3SJed Brown PetscFunctionReturn(0); 622b1a0a8a3SJed Brown } 623b1a0a8a3SJed Brown 624f746d493SDmitry Karpeev static PetscErrorCode PCSetUpOnBlocks_GASM(PC pc) 625b1a0a8a3SJed Brown { 626f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 627b1a0a8a3SJed Brown PetscErrorCode ierr; 628b1a0a8a3SJed Brown PetscInt i; 629b1a0a8a3SJed Brown 630b1a0a8a3SJed Brown PetscFunctionBegin; 631f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 632b1a0a8a3SJed Brown ierr = KSPSetUp(osm->ksp[i]);CHKERRQ(ierr); 633b1a0a8a3SJed Brown } 634b1a0a8a3SJed Brown PetscFunctionReturn(0); 635b1a0a8a3SJed Brown } 636b1a0a8a3SJed Brown 637ea91fabdSFande Kong static PetscErrorCode PCApply_GASM(PC pc,Vec xin,Vec yout) 638b1a0a8a3SJed Brown { 639f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 640b1a0a8a3SJed Brown PetscErrorCode ierr; 641f746d493SDmitry Karpeev PetscInt i; 642ea91fabdSFande Kong Vec x,y; 643b1a0a8a3SJed Brown ScatterMode forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE; 644b1a0a8a3SJed Brown 645b1a0a8a3SJed Brown PetscFunctionBegin; 646ea91fabdSFande Kong if(osm->pctoouter){ 647ea91fabdSFande Kong ierr = VecScatterBegin(osm->pctoouter,xin,osm->pcx,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 648ea91fabdSFande Kong ierr = VecScatterEnd(osm->pctoouter,xin,osm->pcx,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 649ea91fabdSFande Kong x = osm->pcx; 650ea91fabdSFande Kong y = osm->pcy; 651ea91fabdSFande Kong }else{ 652ea91fabdSFande Kong x = xin; 653ea91fabdSFande Kong y = yout; 654ea91fabdSFande Kong } 655b1a0a8a3SJed Brown /* 6566a4f0f73SDmitry Karpeev Support for limiting the restriction or interpolation only to the inner 657b1a0a8a3SJed Brown subdomain values (leaving the other values 0). 658b1a0a8a3SJed Brown */ 659f746d493SDmitry Karpeev if (!(osm->type & PC_GASM_RESTRICT)) { 660b1a0a8a3SJed Brown /* have to zero the work RHS since scatter may leave some slots empty */ 661f746d493SDmitry Karpeev ierr = VecZeroEntries(osm->gx);CHKERRQ(ierr); 6626a4f0f73SDmitry Karpeev ierr = VecScatterBegin(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 6632fa5cd67SKarl Rupp } else { 6646a4f0f73SDmitry Karpeev ierr = VecScatterBegin(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 665b1a0a8a3SJed Brown } 6666a4f0f73SDmitry Karpeev ierr = VecZeroEntries(osm->gy);CHKERRQ(ierr); 6676a4f0f73SDmitry Karpeev if (!(osm->type & PC_GASM_RESTRICT)) { 6686a4f0f73SDmitry Karpeev ierr = VecScatterEnd(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 6692fa5cd67SKarl Rupp } else { 6706a4f0f73SDmitry Karpeev ierr = VecScatterEnd(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 6716a4f0f73SDmitry Karpeev } 67286b96d36SDmitry Karpeev /* do the subdomain solves */ 67386b96d36SDmitry Karpeev for (i=0; i<osm->n; ++i) { 674b1a0a8a3SJed Brown ierr = KSPSolve(osm->ksp[i],osm->x[i],osm->y[i]);CHKERRQ(ierr); 675c0decd05SBarry Smith ierr = KSPCheckSolve(osm->ksp[i],pc,osm->y[i]);CHKERRQ(ierr); 676b1a0a8a3SJed Brown } 677930d09c1SFande Kong /* Do we need to zero y ?? */ 6786a4f0f73SDmitry Karpeev ierr = VecZeroEntries(y);CHKERRQ(ierr); 6796a4f0f73SDmitry Karpeev if (!(osm->type & PC_GASM_INTERPOLATE)) { 6806a4f0f73SDmitry Karpeev ierr = VecScatterBegin(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); 681ea91fabdSFande Kong ierr = VecScatterEnd(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); 6822fa5cd67SKarl Rupp } else { 6836a4f0f73SDmitry Karpeev ierr = VecScatterBegin(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); 684ea91fabdSFande Kong ierr = VecScatterEnd(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); 6856a4f0f73SDmitry Karpeev } 686ea91fabdSFande Kong if(osm->pctoouter){ 687ea91fabdSFande Kong ierr = VecScatterBegin(osm->pctoouter,y,yout,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 688ea91fabdSFande Kong ierr = VecScatterEnd(osm->pctoouter,y,yout,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 689ea91fabdSFande Kong } 690ea91fabdSFande Kong PetscFunctionReturn(0); 691b1a0a8a3SJed Brown } 692b1a0a8a3SJed Brown 693ea91fabdSFande Kong static PetscErrorCode PCApplyTranspose_GASM(PC pc,Vec xin,Vec yout) 694b1a0a8a3SJed Brown { 695f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 696b1a0a8a3SJed Brown PetscErrorCode ierr; 697f746d493SDmitry Karpeev PetscInt i; 698ea91fabdSFande Kong Vec x,y; 699b1a0a8a3SJed Brown ScatterMode forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE; 700b1a0a8a3SJed Brown 701b1a0a8a3SJed Brown PetscFunctionBegin; 702ea91fabdSFande Kong if(osm->pctoouter){ 703ea91fabdSFande Kong ierr = VecScatterBegin(osm->pctoouter,xin,osm->pcx,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 704ea91fabdSFande Kong ierr = VecScatterEnd(osm->pctoouter,xin,osm->pcx,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 705ea91fabdSFande Kong x = osm->pcx; 706ea91fabdSFande Kong y = osm->pcy; 707ea91fabdSFande Kong }else{ 708ea91fabdSFande Kong x = xin; 709ea91fabdSFande Kong y = yout; 710ea91fabdSFande Kong } 711b1a0a8a3SJed Brown /* 712b1a0a8a3SJed Brown Support for limiting the restriction or interpolation to only local 713b1a0a8a3SJed Brown subdomain values (leaving the other values 0). 714b1a0a8a3SJed Brown 715f746d493SDmitry Karpeev Note: these are reversed from the PCApply_GASM() because we are applying the 716b1a0a8a3SJed Brown transpose of the three terms 717b1a0a8a3SJed Brown */ 718f746d493SDmitry Karpeev if (!(osm->type & PC_GASM_INTERPOLATE)) { 719b1a0a8a3SJed Brown /* have to zero the work RHS since scatter may leave some slots empty */ 720f746d493SDmitry Karpeev ierr = VecZeroEntries(osm->gx);CHKERRQ(ierr); 7216a4f0f73SDmitry Karpeev ierr = VecScatterBegin(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 7222fa5cd67SKarl Rupp } else { 7236a4f0f73SDmitry Karpeev ierr = VecScatterBegin(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 724b1a0a8a3SJed Brown } 7256a4f0f73SDmitry Karpeev ierr = VecZeroEntries(osm->gy);CHKERRQ(ierr); 7266a4f0f73SDmitry Karpeev if (!(osm->type & PC_GASM_INTERPOLATE)) { 7276a4f0f73SDmitry Karpeev ierr = VecScatterEnd(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 7282fa5cd67SKarl Rupp } else { 7296a4f0f73SDmitry Karpeev ierr = VecScatterEnd(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 7306a4f0f73SDmitry Karpeev } 731b1a0a8a3SJed Brown /* do the local solves */ 732ab3e923fSDmitry 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. */ 733b1a0a8a3SJed Brown ierr = KSPSolveTranspose(osm->ksp[i],osm->x[i],osm->y[i]);CHKERRQ(ierr); 734c0decd05SBarry Smith ierr = KSPCheckSolve(osm->ksp[i],pc,osm->y[i]);CHKERRQ(ierr); 735b1a0a8a3SJed Brown } 7366a4f0f73SDmitry Karpeev ierr = VecZeroEntries(y);CHKERRQ(ierr); 7376a4f0f73SDmitry Karpeev if (!(osm->type & PC_GASM_RESTRICT)) { 7386a4f0f73SDmitry Karpeev ierr = VecScatterBegin(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); 7396a4f0f73SDmitry Karpeev ierr = VecScatterEnd(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); 7402fa5cd67SKarl Rupp } else { 7416a4f0f73SDmitry Karpeev ierr = VecScatterBegin(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); 7426a4f0f73SDmitry Karpeev ierr = VecScatterEnd(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); 7436a4f0f73SDmitry Karpeev } 744ea91fabdSFande Kong if(osm->pctoouter){ 745ea91fabdSFande Kong ierr = VecScatterBegin(osm->pctoouter,y,yout,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 746ea91fabdSFande Kong ierr = VecScatterEnd(osm->pctoouter,y,yout,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 747ea91fabdSFande Kong } 748b1a0a8a3SJed Brown PetscFunctionReturn(0); 749b1a0a8a3SJed Brown } 750b1a0a8a3SJed Brown 751a06653b4SBarry Smith static PetscErrorCode PCReset_GASM(PC pc) 752b1a0a8a3SJed Brown { 753f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 754b1a0a8a3SJed Brown PetscErrorCode ierr; 755b1a0a8a3SJed Brown PetscInt i; 756b1a0a8a3SJed Brown 757b1a0a8a3SJed Brown PetscFunctionBegin; 758b1a0a8a3SJed Brown if (osm->ksp) { 759f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 760a06653b4SBarry Smith ierr = KSPReset(osm->ksp[i]);CHKERRQ(ierr); 761b1a0a8a3SJed Brown } 762b1a0a8a3SJed Brown } 763b1a0a8a3SJed Brown if (osm->pmat) { 764f746d493SDmitry Karpeev if (osm->n > 0) { 765df750dc8SHong Zhang PetscMPIInt size; 766df750dc8SHong Zhang ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr); 767df750dc8SHong Zhang if (size > 1) { 7687dae84e0SHong Zhang /* osm->pmat is created by MatCreateSubMatricesMPI(), cannot use MatDestroySubMatrices() */ 769ab3e923fSDmitry Karpeev ierr = MatDestroyMatrices(osm->n,&osm->pmat);CHKERRQ(ierr); 770df750dc8SHong Zhang } else { 771df750dc8SHong Zhang ierr = MatDestroySubMatrices(osm->n,&osm->pmat);CHKERRQ(ierr); 772df750dc8SHong Zhang } 773b1a0a8a3SJed Brown } 774b1a0a8a3SJed Brown } 775d34fcf5fSBarry Smith if (osm->x) { 776f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 777fcfd50ebSBarry Smith ierr = VecDestroy(&osm->x[i]);CHKERRQ(ierr); 778fcfd50ebSBarry Smith ierr = VecDestroy(&osm->y[i]);CHKERRQ(ierr); 779b1a0a8a3SJed Brown } 780d34fcf5fSBarry Smith } 7816bf464f9SBarry Smith ierr = VecDestroy(&osm->gx);CHKERRQ(ierr); 7826bf464f9SBarry Smith ierr = VecDestroy(&osm->gy);CHKERRQ(ierr); 783ab3e923fSDmitry Karpeev 7846a4f0f73SDmitry Karpeev ierr = VecScatterDestroy(&osm->gorestriction);CHKERRQ(ierr); 7856a4f0f73SDmitry Karpeev ierr = VecScatterDestroy(&osm->girestriction);CHKERRQ(ierr); 7868f3b4b4dSDmitry Karpeev if (!osm->user_subdomains) { 7872c112581SDmitry Karpeev ierr = PCGASMDestroySubdomains(osm->n,&osm->ois,&osm->iis);CHKERRQ(ierr); 7888f3b4b4dSDmitry Karpeev osm->N = PETSC_DETERMINE; 7898f3b4b4dSDmitry Karpeev osm->nmax = PETSC_DETERMINE; 7908f3b4b4dSDmitry Karpeev } 791ea91fabdSFande Kong if(osm->pctoouter){ 792ea91fabdSFande Kong ierr = VecScatterDestroy(&(osm->pctoouter));CHKERRQ(ierr); 793ea91fabdSFande Kong } 794ea91fabdSFande Kong if(osm->permutationIS){ 795ea91fabdSFande Kong ierr = ISDestroy(&(osm->permutationIS));CHKERRQ(ierr); 796ea91fabdSFande Kong } 797ea91fabdSFande Kong if(osm->pcx){ 798ea91fabdSFande Kong ierr = VecDestroy(&(osm->pcx));CHKERRQ(ierr); 799ea91fabdSFande Kong } 800ea91fabdSFande Kong if(osm->pcy){ 801ea91fabdSFande Kong ierr = VecDestroy(&(osm->pcy));CHKERRQ(ierr); 802ea91fabdSFande Kong } 803ea91fabdSFande Kong if(osm->permutationP){ 804ea91fabdSFande Kong ierr = MatDestroy(&(osm->permutationP));CHKERRQ(ierr); 805ea91fabdSFande Kong } 806ea91fabdSFande Kong if(osm->pcmat){ 807ea91fabdSFande Kong ierr = MatDestroy(&osm->pcmat);CHKERRQ(ierr); 808ea91fabdSFande Kong } 809a06653b4SBarry Smith PetscFunctionReturn(0); 810a06653b4SBarry Smith } 811a06653b4SBarry Smith 812a06653b4SBarry Smith static PetscErrorCode PCDestroy_GASM(PC pc) 813a06653b4SBarry Smith { 814a06653b4SBarry Smith PC_GASM *osm = (PC_GASM*)pc->data; 815a06653b4SBarry Smith PetscErrorCode ierr; 816a06653b4SBarry Smith PetscInt i; 817a06653b4SBarry Smith 818a06653b4SBarry Smith PetscFunctionBegin; 819a06653b4SBarry Smith ierr = PCReset_GASM(pc);CHKERRQ(ierr); 820135757f6SDmitry Karpeev /* PCReset will not destroy subdomains, if user_subdomains is true. */ 821135757f6SDmitry Karpeev ierr = PCGASMDestroySubdomains(osm->n,&osm->ois,&osm->iis);CHKERRQ(ierr); 822a06653b4SBarry Smith if (osm->ksp) { 823a06653b4SBarry Smith for (i=0; i<osm->n; i++) { 8246bf464f9SBarry Smith ierr = KSPDestroy(&osm->ksp[i]);CHKERRQ(ierr); 825a06653b4SBarry Smith } 826a06653b4SBarry Smith ierr = PetscFree(osm->ksp);CHKERRQ(ierr); 827a06653b4SBarry Smith } 828a06653b4SBarry Smith ierr = PetscFree(osm->x);CHKERRQ(ierr); 829a06653b4SBarry Smith ierr = PetscFree(osm->y);CHKERRQ(ierr); 830c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 831b1a0a8a3SJed Brown PetscFunctionReturn(0); 832b1a0a8a3SJed Brown } 833b1a0a8a3SJed Brown 8344416b707SBarry Smith static PetscErrorCode PCSetFromOptions_GASM(PetscOptionItems *PetscOptionsObject,PC pc) 835a6dfd86eSKarl Rupp { 836f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 837b1a0a8a3SJed Brown PetscErrorCode ierr; 838b1a0a8a3SJed Brown PetscInt blocks,ovl; 83957501b6eSBarry Smith PetscBool flg; 840f746d493SDmitry Karpeev PCGASMType gasmtype; 841b1a0a8a3SJed Brown 842b1a0a8a3SJed Brown PetscFunctionBegin; 843e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"Generalized additive Schwarz options");CHKERRQ(ierr); 8448f3b4b4dSDmitry Karpeev ierr = PetscOptionsBool("-pc_gasm_use_dm_subdomains","If subdomains aren't set, use DMCreateDomainDecomposition() to define subdomains.","PCGASMSetUseDMSubdomains",osm->dm_subdomains,&osm->dm_subdomains,&flg);CHKERRQ(ierr); 8458f3b4b4dSDmitry Karpeev ierr = PetscOptionsInt("-pc_gasm_total_subdomains","Total number of subdomains across communicator","PCGASMSetTotalSubdomains",osm->N,&blocks,&flg);CHKERRQ(ierr); 84665db9045SDmitry Karpeev if (flg) { 8478f3b4b4dSDmitry Karpeev ierr = PCGASMSetTotalSubdomains(pc,blocks);CHKERRQ(ierr); 84865db9045SDmitry Karpeev } 84906b43e7eSDmitry Karpeev ierr = PetscOptionsInt("-pc_gasm_overlap","Number of overlapping degrees of freedom","PCGASMSetOverlap",osm->overlap,&ovl,&flg);CHKERRQ(ierr); 85065db9045SDmitry Karpeev if (flg) { 85165db9045SDmitry Karpeev ierr = PCGASMSetOverlap(pc,ovl);CHKERRQ(ierr); 852d709ea83SDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 85365db9045SDmitry Karpeev } 854b1a0a8a3SJed Brown flg = PETSC_FALSE; 855f746d493SDmitry Karpeev ierr = PetscOptionsEnum("-pc_gasm_type","Type of restriction/extension","PCGASMSetType",PCGASMTypes,(PetscEnum)osm->type,(PetscEnum*)&gasmtype,&flg);CHKERRQ(ierr); 856f746d493SDmitry Karpeev if (flg) {ierr = PCGASMSetType(pc,gasmtype);CHKERRQ(ierr);} 857ea91fabdSFande Kong ierr = PetscOptionsBool("-pc_gasm_use_hierachical_partitioning","use hierarchical partitioning",NULL,osm->hierarchicalpartitioning,&osm->hierarchicalpartitioning,&flg);CHKERRQ(ierr); 858b1a0a8a3SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 859b1a0a8a3SJed Brown PetscFunctionReturn(0); 860b1a0a8a3SJed Brown } 861b1a0a8a3SJed Brown 862b1a0a8a3SJed Brown /*------------------------------------------------------------------------------------*/ 863b1a0a8a3SJed Brown 8648f3b4b4dSDmitry Karpeev /*@ 8658f3b4b4dSDmitry Karpeev PCGASMSetTotalSubdomains - sets the total number of subdomains to use across the 8668f3b4b4dSDmitry Karpeev communicator. 8678f3b4b4dSDmitry Karpeev Logically collective on pc 8688f3b4b4dSDmitry Karpeev 8698f3b4b4dSDmitry Karpeev Input Parameters: 8708f3b4b4dSDmitry Karpeev + pc - the preconditioner 8718f3b4b4dSDmitry Karpeev - N - total number of subdomains 8728f3b4b4dSDmitry Karpeev 8738f3b4b4dSDmitry Karpeev 8748f3b4b4dSDmitry Karpeev Level: beginner 8758f3b4b4dSDmitry Karpeev 8768f3b4b4dSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMSetOverlap() 8778f3b4b4dSDmitry Karpeev PCGASMCreateSubdomains2D() 8788f3b4b4dSDmitry Karpeev @*/ 8798f3b4b4dSDmitry Karpeev PetscErrorCode PCGASMSetTotalSubdomains(PC pc,PetscInt N) 8808f3b4b4dSDmitry Karpeev { 8818f3b4b4dSDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 8828f3b4b4dSDmitry Karpeev PetscMPIInt size,rank; 8838f3b4b4dSDmitry Karpeev PetscErrorCode ierr; 8848f3b4b4dSDmitry Karpeev 8858f3b4b4dSDmitry Karpeev PetscFunctionBegin; 8868f3b4b4dSDmitry Karpeev if (N < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Total number of subdomains must be 1 or more, got N = %D",N); 8878f3b4b4dSDmitry Karpeev if (pc->setupcalled) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetTotalSubdomains() should be called before calling PCSetUp()."); 8888f3b4b4dSDmitry Karpeev 8892c112581SDmitry Karpeev ierr = PCGASMDestroySubdomains(osm->n,&osm->iis,&osm->ois);CHKERRQ(ierr); 8908f3b4b4dSDmitry Karpeev osm->ois = osm->iis = NULL; 8918f3b4b4dSDmitry Karpeev 8928f3b4b4dSDmitry Karpeev ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr); 8938f3b4b4dSDmitry Karpeev ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);CHKERRQ(ierr); 8948f3b4b4dSDmitry Karpeev osm->N = N; 8958f3b4b4dSDmitry Karpeev osm->n = PETSC_DETERMINE; 8968f3b4b4dSDmitry Karpeev osm->nmax = PETSC_DETERMINE; 8978f3b4b4dSDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 8988f3b4b4dSDmitry Karpeev PetscFunctionReturn(0); 8998f3b4b4dSDmitry Karpeev } 9008f3b4b4dSDmitry Karpeev 901670417b2SFande Kong 902b541e6a4SDmitry Karpeev static PetscErrorCode PCGASMSetSubdomains_GASM(PC pc,PetscInt n,IS iis[],IS ois[]) 903b1a0a8a3SJed Brown { 904f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 905b1a0a8a3SJed Brown PetscErrorCode ierr; 906b1a0a8a3SJed Brown PetscInt i; 907b1a0a8a3SJed Brown 908b1a0a8a3SJed Brown PetscFunctionBegin; 9098f3b4b4dSDmitry Karpeev if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Each process must have 1 or more subdomains, got n = %D",n); 9108f3b4b4dSDmitry Karpeev if (pc->setupcalled) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetSubdomains() should be called before calling PCSetUp()."); 911b1a0a8a3SJed Brown 9122c112581SDmitry Karpeev ierr = PCGASMDestroySubdomains(osm->n,&osm->iis,&osm->ois);CHKERRQ(ierr); 9138f3b4b4dSDmitry Karpeev osm->iis = osm->ois = NULL; 9148f3b4b4dSDmitry Karpeev osm->n = n; 9158f3b4b4dSDmitry Karpeev osm->N = PETSC_DETERMINE; 9168f3b4b4dSDmitry Karpeev osm->nmax = PETSC_DETERMINE; 917a35b7b57SDmitry Karpeev if (ois) { 918785e854fSJed Brown ierr = PetscMalloc1(n,&osm->ois);CHKERRQ(ierr); 9198f3b4b4dSDmitry Karpeev for (i=0; i<n; i++) { 9208f3b4b4dSDmitry Karpeev ierr = PetscObjectReference((PetscObject)ois[i]);CHKERRQ(ierr); 9218f3b4b4dSDmitry Karpeev osm->ois[i] = ois[i]; 9228f3b4b4dSDmitry Karpeev } 9238f3b4b4dSDmitry Karpeev /* 9248f3b4b4dSDmitry Karpeev Since the user set the outer subdomains, even if nontrivial overlap was requested via PCGASMSetOverlap(), 9258f3b4b4dSDmitry Karpeev it will be ignored. To avoid confusion later on (e.g., when viewing the PC), the overlap size is set to -1. 9268f3b4b4dSDmitry Karpeev */ 927b1a0a8a3SJed Brown osm->overlap = -1; 928670417b2SFande Kong /* inner subdomains must be provided */ 929670417b2SFande Kong if (!iis) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"inner indices have to be provided \n"); 930670417b2SFande Kong }/* end if */ 931a35b7b57SDmitry Karpeev if (iis) { 932785e854fSJed Brown ierr = PetscMalloc1(n,&osm->iis);CHKERRQ(ierr); 9338f3b4b4dSDmitry Karpeev for (i=0; i<n; i++) { 9348f3b4b4dSDmitry Karpeev ierr = PetscObjectReference((PetscObject)iis[i]);CHKERRQ(ierr); 9358f3b4b4dSDmitry Karpeev osm->iis[i] = iis[i]; 9368f3b4b4dSDmitry Karpeev } 937a35b7b57SDmitry Karpeev if (!ois) { 938390e1bf2SBarry Smith osm->ois = NULL; 939670417b2SFande Kong /* if user does not provide outer indices, we will create the corresponding outer indices using osm->overlap =1 in PCSetUp_GASM */ 940670417b2SFande Kong } 941670417b2SFande Kong } 94276bd3646SJed Brown if (PetscDefined(USE_DEBUG)) { 943670417b2SFande Kong PetscInt j,rstart,rend,*covered,lsize; 944670417b2SFande Kong const PetscInt *indices; 945670417b2SFande Kong /* check if the inner indices cover and only cover the local portion of the preconditioning matrix */ 946670417b2SFande Kong ierr = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr); 947670417b2SFande Kong ierr = PetscCalloc1(rend-rstart,&covered);CHKERRQ(ierr); 948670417b2SFande Kong /* check if the current processor owns indices from others */ 949a35b7b57SDmitry Karpeev for (i=0; i<n; i++) { 950670417b2SFande Kong ierr = ISGetIndices(osm->iis[i],&indices);CHKERRQ(ierr); 951670417b2SFande Kong ierr = ISGetLocalSize(osm->iis[i],&lsize);CHKERRQ(ierr); 952670417b2SFande Kong for (j=0; j<lsize; j++) { 953670417b2SFande Kong if (indices[j]<rstart || indices[j]>=rend) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"inner subdomains can not own an index %d from other processors", indices[j]); 954670417b2SFande Kong else if (covered[indices[j]-rstart]==1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"inner subdomains can not have an overlapping index %d ",indices[j]); 955670417b2SFande Kong else covered[indices[j]-rstart] = 1; 956a35b7b57SDmitry Karpeev } 957670417b2SFande Kong ierr = ISRestoreIndices(osm->iis[i],&indices);CHKERRQ(ierr); 9588f3b4b4dSDmitry Karpeev } 959670417b2SFande Kong /* check if we miss any indices */ 960670417b2SFande Kong for (i=rstart; i<rend; i++) { 961670417b2SFande Kong if (!covered[i-rstart]) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"local entity %d was not covered by inner subdomains",i); 962a35b7b57SDmitry Karpeev } 963670417b2SFande Kong ierr = PetscFree(covered);CHKERRQ(ierr); 964a35b7b57SDmitry Karpeev } 9658f3b4b4dSDmitry Karpeev if (iis) osm->user_subdomains = PETSC_TRUE; 966b1a0a8a3SJed Brown PetscFunctionReturn(0); 967b1a0a8a3SJed Brown } 968b1a0a8a3SJed Brown 969b1a0a8a3SJed Brown 970f7a08781SBarry Smith static PetscErrorCode PCGASMSetOverlap_GASM(PC pc,PetscInt ovl) 971b1a0a8a3SJed Brown { 972f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 973b1a0a8a3SJed Brown 974b1a0a8a3SJed Brown PetscFunctionBegin; 975ce94432eSBarry Smith if (ovl < 0) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap value requested"); 976ce94432eSBarry Smith if (pc->setupcalled && ovl != osm->overlap) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetOverlap() should be called before PCSetUp()."); 9772fa5cd67SKarl Rupp if (!pc->setupcalled) osm->overlap = ovl; 978b1a0a8a3SJed Brown PetscFunctionReturn(0); 979b1a0a8a3SJed Brown } 980b1a0a8a3SJed Brown 981f7a08781SBarry Smith static PetscErrorCode PCGASMSetType_GASM(PC pc,PCGASMType type) 982b1a0a8a3SJed Brown { 983f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 984b1a0a8a3SJed Brown 985b1a0a8a3SJed Brown PetscFunctionBegin; 986b1a0a8a3SJed Brown osm->type = type; 987b1a0a8a3SJed Brown osm->type_set = PETSC_TRUE; 988b1a0a8a3SJed Brown PetscFunctionReturn(0); 989b1a0a8a3SJed Brown } 990b1a0a8a3SJed Brown 991f7a08781SBarry Smith static PetscErrorCode PCGASMSetSortIndices_GASM(PC pc,PetscBool doSort) 992b1a0a8a3SJed Brown { 993f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 994b1a0a8a3SJed Brown 995b1a0a8a3SJed Brown PetscFunctionBegin; 996b1a0a8a3SJed Brown osm->sort_indices = doSort; 997b1a0a8a3SJed Brown PetscFunctionReturn(0); 998b1a0a8a3SJed Brown } 999b1a0a8a3SJed Brown 100044fe18e5SDmitry Karpeev /* 10018f3b4b4dSDmitry Karpeev FIXME: This routine might need to be modified now that multiple ranks per subdomain are allowed. 100244fe18e5SDmitry Karpeev In particular, it would upset the global subdomain number calculation. 100344fe18e5SDmitry Karpeev */ 1004f7a08781SBarry Smith static PetscErrorCode PCGASMGetSubKSP_GASM(PC pc,PetscInt *n,PetscInt *first,KSP **ksp) 1005b1a0a8a3SJed Brown { 1006f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 1007b1a0a8a3SJed Brown PetscErrorCode ierr; 1008b1a0a8a3SJed Brown 1009b1a0a8a3SJed Brown PetscFunctionBegin; 101034a84908Sprj- 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"); 1011b1a0a8a3SJed Brown 10122fa5cd67SKarl Rupp if (n) *n = osm->n; 1013ab3e923fSDmitry Karpeev if (first) { 1014ce94432eSBarry Smith ierr = MPI_Scan(&osm->n,first,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr); 1015ab3e923fSDmitry Karpeev *first -= osm->n; 1016b1a0a8a3SJed Brown } 1017b1a0a8a3SJed Brown if (ksp) { 1018b1a0a8a3SJed Brown /* Assume that local solves are now different; not necessarily 101906b43e7eSDmitry Karpeev true, though! This flag is used only for PCView_GASM() */ 1020b1a0a8a3SJed Brown *ksp = osm->ksp; 10216a4f0f73SDmitry Karpeev osm->same_subdomain_solvers = PETSC_FALSE; 1022b1a0a8a3SJed Brown } 1023b1a0a8a3SJed Brown PetscFunctionReturn(0); 1024ab3e923fSDmitry Karpeev } /* PCGASMGetSubKSP_GASM() */ 1025b1a0a8a3SJed Brown 1026b1a0a8a3SJed Brown /*@C 102706b43e7eSDmitry Karpeev PCGASMSetSubdomains - Sets the subdomains for this processor 102806b43e7eSDmitry Karpeev for the additive Schwarz preconditioner. 1029b1a0a8a3SJed Brown 10308f3b4b4dSDmitry Karpeev Collective on pc 1031b1a0a8a3SJed Brown 1032b1a0a8a3SJed Brown Input Parameters: 10338f3b4b4dSDmitry Karpeev + pc - the preconditioner object 103406b43e7eSDmitry Karpeev . n - the number of subdomains for this processor 10358f3b4b4dSDmitry Karpeev . iis - the index sets that define the inner subdomains (or NULL for PETSc to determine subdomains) 10368f3b4b4dSDmitry Karpeev - ois - the index sets that define the outer subdomains (or NULL to use the same as iis, or to construct by expanding iis by the requested overlap) 1037b1a0a8a3SJed Brown 1038b1a0a8a3SJed Brown Notes: 103906b43e7eSDmitry Karpeev The IS indices use the parallel, global numbering of the vector entries. 10406a4f0f73SDmitry Karpeev Inner subdomains are those where the correction is applied. 10416a4f0f73SDmitry Karpeev Outer subdomains are those where the residual necessary to obtain the 10426a4f0f73SDmitry Karpeev corrections is obtained (see PCGASMType for the use of inner/outer subdomains). 10436a4f0f73SDmitry Karpeev Both inner and outer subdomains can extend over several processors. 10446a4f0f73SDmitry Karpeev This processor's portion of a subdomain is known as a local subdomain. 10456a4f0f73SDmitry Karpeev 1046670417b2SFande Kong Inner subdomains can not overlap with each other, do not have any entities from remote processors, 1047670417b2SFande Kong and have to cover the entire local subdomain owned by the current processor. The index sets on each 1048670417b2SFande Kong process should be ordered such that the ith local subdomain is connected to the ith remote subdomain 1049670417b2SFande Kong on another MPI process. 1050670417b2SFande Kong 10516a4f0f73SDmitry Karpeev By default the GASM preconditioner uses 1 (local) subdomain per processor. 10526a4f0f73SDmitry Karpeev 1053b1a0a8a3SJed Brown 1054b1a0a8a3SJed Brown Level: advanced 1055b1a0a8a3SJed Brown 10568f3b4b4dSDmitry Karpeev .seealso: PCGASMSetNumSubdomains(), PCGASMSetOverlap(), PCGASMGetSubKSP(), 105706b43e7eSDmitry Karpeev PCGASMCreateSubdomains2D(), PCGASMGetSubdomains() 1058b1a0a8a3SJed Brown @*/ 10596a4f0f73SDmitry Karpeev PetscErrorCode PCGASMSetSubdomains(PC pc,PetscInt n,IS iis[],IS ois[]) 1060b1a0a8a3SJed Brown { 10618f3b4b4dSDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 1062b1a0a8a3SJed Brown PetscErrorCode ierr; 1063b1a0a8a3SJed Brown 1064b1a0a8a3SJed Brown PetscFunctionBegin; 1065b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 10666a4f0f73SDmitry Karpeev ierr = PetscTryMethod(pc,"PCGASMSetSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,n,iis,ois));CHKERRQ(ierr); 10678f3b4b4dSDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 1068b1a0a8a3SJed Brown PetscFunctionReturn(0); 1069b1a0a8a3SJed Brown } 1070b1a0a8a3SJed Brown 1071b1a0a8a3SJed Brown 1072b1a0a8a3SJed Brown /*@ 1073f746d493SDmitry Karpeev PCGASMSetOverlap - Sets the overlap between a pair of subdomains for the 1074b1a0a8a3SJed Brown additive Schwarz preconditioner. Either all or no processors in the 10758f3b4b4dSDmitry Karpeev pc communicator must call this routine. 1076b1a0a8a3SJed Brown 10778f3b4b4dSDmitry Karpeev Logically Collective on pc 1078b1a0a8a3SJed Brown 1079b1a0a8a3SJed Brown Input Parameters: 1080b1a0a8a3SJed Brown + pc - the preconditioner context 10818f3b4b4dSDmitry Karpeev - ovl - the amount of overlap between subdomains (ovl >= 0, default value = 0) 1082b1a0a8a3SJed Brown 1083b1a0a8a3SJed Brown Options Database Key: 108406b43e7eSDmitry Karpeev . -pc_gasm_overlap <overlap> - Sets overlap 1085b1a0a8a3SJed Brown 1086b1a0a8a3SJed Brown Notes: 108706b43e7eSDmitry Karpeev By default the GASM preconditioner uses 1 subdomain per processor. To use 10888f3b4b4dSDmitry Karpeev multiple subdomain per perocessor or "straddling" subdomains that intersect 10898f3b4b4dSDmitry Karpeev multiple processors use PCGASMSetSubdomains() (or option -pc_gasm_total_subdomains <n>). 1090b1a0a8a3SJed Brown 10918f3b4b4dSDmitry Karpeev The overlap defaults to 0, so if one desires that no additional 1092b1a0a8a3SJed Brown overlap be computed beyond what may have been set with a call to 10938f3b4b4dSDmitry Karpeev PCGASMSetSubdomains(), then ovl must be set to be 0. In particular, if one does 10948f3b4b4dSDmitry Karpeev not explicitly set the subdomains in application code, then all overlap would be computed 1095f746d493SDmitry Karpeev internally by PETSc, and using an overlap of 0 would result in an GASM 1096b1a0a8a3SJed Brown variant that is equivalent to the block Jacobi preconditioner. 1097b1a0a8a3SJed Brown 1098b1a0a8a3SJed Brown Note that one can define initial index sets with any overlap via 109906b43e7eSDmitry Karpeev PCGASMSetSubdomains(); the routine PCGASMSetOverlap() merely allows 110006b43e7eSDmitry Karpeev PETSc to extend that overlap further, if desired. 1101b1a0a8a3SJed Brown 1102b1a0a8a3SJed Brown Level: intermediate 1103b1a0a8a3SJed Brown 11048f3b4b4dSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMGetSubKSP(), 110506b43e7eSDmitry Karpeev PCGASMCreateSubdomains2D(), PCGASMGetSubdomains() 1106b1a0a8a3SJed Brown @*/ 11077087cfbeSBarry Smith PetscErrorCode PCGASMSetOverlap(PC pc,PetscInt ovl) 1108b1a0a8a3SJed Brown { 1109b1a0a8a3SJed Brown PetscErrorCode ierr; 11108f3b4b4dSDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 1111b1a0a8a3SJed Brown 1112b1a0a8a3SJed Brown PetscFunctionBegin; 1113b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1114b1a0a8a3SJed Brown PetscValidLogicalCollectiveInt(pc,ovl,2); 1115f746d493SDmitry Karpeev ierr = PetscTryMethod(pc,"PCGASMSetOverlap_C",(PC,PetscInt),(pc,ovl));CHKERRQ(ierr); 11168f3b4b4dSDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 1117b1a0a8a3SJed Brown PetscFunctionReturn(0); 1118b1a0a8a3SJed Brown } 1119b1a0a8a3SJed Brown 1120b1a0a8a3SJed Brown /*@ 1121f746d493SDmitry Karpeev PCGASMSetType - Sets the type of restriction and interpolation used 1122b1a0a8a3SJed Brown for local problems in the additive Schwarz method. 1123b1a0a8a3SJed Brown 1124b1a0a8a3SJed Brown Logically Collective on PC 1125b1a0a8a3SJed Brown 1126b1a0a8a3SJed Brown Input Parameters: 1127b1a0a8a3SJed Brown + pc - the preconditioner context 1128f746d493SDmitry Karpeev - type - variant of GASM, one of 1129b1a0a8a3SJed Brown .vb 1130f746d493SDmitry Karpeev PC_GASM_BASIC - full interpolation and restriction 1131f746d493SDmitry Karpeev PC_GASM_RESTRICT - full restriction, local processor interpolation 1132f746d493SDmitry Karpeev PC_GASM_INTERPOLATE - full interpolation, local processor restriction 1133f746d493SDmitry Karpeev PC_GASM_NONE - local processor restriction and interpolation 1134b1a0a8a3SJed Brown .ve 1135b1a0a8a3SJed Brown 1136b1a0a8a3SJed Brown Options Database Key: 1137f746d493SDmitry Karpeev . -pc_gasm_type [basic,restrict,interpolate,none] - Sets GASM type 1138b1a0a8a3SJed Brown 1139b1a0a8a3SJed Brown Level: intermediate 1140b1a0a8a3SJed Brown 11418f3b4b4dSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMGetSubKSP(), 1142f746d493SDmitry Karpeev PCGASMCreateSubdomains2D() 1143b1a0a8a3SJed Brown @*/ 11447087cfbeSBarry Smith PetscErrorCode PCGASMSetType(PC pc,PCGASMType type) 1145b1a0a8a3SJed Brown { 1146b1a0a8a3SJed Brown PetscErrorCode ierr; 1147b1a0a8a3SJed Brown 1148b1a0a8a3SJed Brown PetscFunctionBegin; 1149b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1150b1a0a8a3SJed Brown PetscValidLogicalCollectiveEnum(pc,type,2); 1151f746d493SDmitry Karpeev ierr = PetscTryMethod(pc,"PCGASMSetType_C",(PC,PCGASMType),(pc,type));CHKERRQ(ierr); 1152b1a0a8a3SJed Brown PetscFunctionReturn(0); 1153b1a0a8a3SJed Brown } 1154b1a0a8a3SJed Brown 1155b1a0a8a3SJed Brown /*@ 1156f746d493SDmitry Karpeev PCGASMSetSortIndices - Determines whether subdomain indices are sorted. 1157b1a0a8a3SJed Brown 1158b1a0a8a3SJed Brown Logically Collective on PC 1159b1a0a8a3SJed Brown 1160b1a0a8a3SJed Brown Input Parameters: 1161b1a0a8a3SJed Brown + pc - the preconditioner context 1162b1a0a8a3SJed Brown - doSort - sort the subdomain indices 1163b1a0a8a3SJed Brown 1164b1a0a8a3SJed Brown Level: intermediate 1165b1a0a8a3SJed Brown 11668f3b4b4dSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMGetSubKSP(), 1167f746d493SDmitry Karpeev PCGASMCreateSubdomains2D() 1168b1a0a8a3SJed Brown @*/ 11697087cfbeSBarry Smith PetscErrorCode PCGASMSetSortIndices(PC pc,PetscBool doSort) 1170b1a0a8a3SJed Brown { 1171b1a0a8a3SJed Brown PetscErrorCode ierr; 1172b1a0a8a3SJed Brown 1173b1a0a8a3SJed Brown PetscFunctionBegin; 1174b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1175b1a0a8a3SJed Brown PetscValidLogicalCollectiveBool(pc,doSort,2); 1176f746d493SDmitry Karpeev ierr = PetscTryMethod(pc,"PCGASMSetSortIndices_C",(PC,PetscBool),(pc,doSort));CHKERRQ(ierr); 1177b1a0a8a3SJed Brown PetscFunctionReturn(0); 1178b1a0a8a3SJed Brown } 1179b1a0a8a3SJed Brown 1180b1a0a8a3SJed Brown /*@C 1181f746d493SDmitry Karpeev PCGASMGetSubKSP - Gets the local KSP contexts for all blocks on 1182b1a0a8a3SJed Brown this processor. 1183b1a0a8a3SJed Brown 1184b1a0a8a3SJed Brown Collective on PC iff first_local is requested 1185b1a0a8a3SJed Brown 1186b1a0a8a3SJed Brown Input Parameter: 1187b1a0a8a3SJed Brown . pc - the preconditioner context 1188b1a0a8a3SJed Brown 1189b1a0a8a3SJed Brown Output Parameters: 11900298fd71SBarry Smith + n_local - the number of blocks on this processor or NULL 11910298fd71SBarry Smith . first_local - the global number of the first block on this processor or NULL, 11920298fd71SBarry Smith all processors must request or all must pass NULL 1193b1a0a8a3SJed Brown - ksp - the array of KSP contexts 1194b1a0a8a3SJed Brown 1195b1a0a8a3SJed Brown Note: 1196f746d493SDmitry Karpeev After PCGASMGetSubKSP() the array of KSPes is not to be freed 1197b1a0a8a3SJed Brown 1198b1a0a8a3SJed Brown Currently for some matrix implementations only 1 block per processor 1199b1a0a8a3SJed Brown is supported. 1200b1a0a8a3SJed Brown 1201f746d493SDmitry Karpeev You must call KSPSetUp() before calling PCGASMGetSubKSP(). 1202b1a0a8a3SJed Brown 1203b1a0a8a3SJed Brown Level: advanced 1204b1a0a8a3SJed Brown 12058f3b4b4dSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMSetOverlap(), 1206f746d493SDmitry Karpeev PCGASMCreateSubdomains2D(), 1207b1a0a8a3SJed Brown @*/ 12087087cfbeSBarry Smith PetscErrorCode PCGASMGetSubKSP(PC pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp[]) 1209b1a0a8a3SJed Brown { 1210b1a0a8a3SJed Brown PetscErrorCode ierr; 1211b1a0a8a3SJed Brown 1212b1a0a8a3SJed Brown PetscFunctionBegin; 1213b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1214f746d493SDmitry Karpeev ierr = PetscUseMethod(pc,"PCGASMGetSubKSP_C",(PC,PetscInt*,PetscInt*,KSP **),(pc,n_local,first_local,ksp));CHKERRQ(ierr); 1215b1a0a8a3SJed Brown PetscFunctionReturn(0); 1216b1a0a8a3SJed Brown } 1217b1a0a8a3SJed Brown 1218b1a0a8a3SJed Brown /* -------------------------------------------------------------------------------------*/ 1219b1a0a8a3SJed Brown /*MC 1220f746d493SDmitry Karpeev PCGASM - Use the (restricted) additive Schwarz method, each block is (approximately) solved with 1221b1a0a8a3SJed Brown its own KSP object. 1222b1a0a8a3SJed Brown 1223b1a0a8a3SJed Brown Options Database Keys: 12248f3b4b4dSDmitry Karpeev + -pc_gasm_total_subdomains <n> - Sets total number of local subdomains to be distributed among processors 122506b43e7eSDmitry Karpeev . -pc_gasm_view_subdomains - activates the printing of subdomain indices in PCView(), -ksp_view or -snes_view 122606b43e7eSDmitry Karpeev . -pc_gasm_print_subdomains - activates the printing of subdomain indices in PCSetUp() 122706b43e7eSDmitry Karpeev . -pc_gasm_overlap <ovl> - Sets overlap by which to (automatically) extend local subdomains 1228f746d493SDmitry Karpeev - -pc_gasm_type [basic,restrict,interpolate,none] - Sets GASM type 1229b1a0a8a3SJed Brown 1230b1a0a8a3SJed Brown IMPORTANT: If you run with, for example, 3 blocks on 1 processor or 3 blocks on 3 processors you 1231f746d493SDmitry Karpeev will get a different convergence rate due to the default option of -pc_gasm_type restrict. Use 1232f746d493SDmitry Karpeev -pc_gasm_type basic to use the standard GASM. 1233b1a0a8a3SJed Brown 123495452b02SPatrick Sanan Notes: 123595452b02SPatrick Sanan Blocks can be shared by multiple processes. 1236b1a0a8a3SJed Brown 1237b1a0a8a3SJed Brown To set options on the solvers for each block append -sub_ to all the KSP, and PC 1238b1a0a8a3SJed Brown options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly 1239b1a0a8a3SJed Brown 1240f746d493SDmitry Karpeev To set the options on the solvers separate for each block call PCGASMGetSubKSP() 1241b1a0a8a3SJed Brown and set the options directly on the resulting KSP object (you can access its PC 1242b1a0a8a3SJed Brown with KSPGetPC()) 1243b1a0a8a3SJed Brown 1244b1a0a8a3SJed Brown 1245b1a0a8a3SJed Brown Level: beginner 1246b1a0a8a3SJed Brown 1247b1a0a8a3SJed Brown References: 124896a0c994SBarry Smith + 1. - M Dryja, OB Widlund, An additive variant of the Schwarz alternating method for the case of many subregions 124996a0c994SBarry Smith Courant Institute, New York University Technical report 125096a0c994SBarry Smith - 2. - Barry Smith, Petter Bjorstad, and William Gropp, Domain Decompositions: Parallel Multilevel Methods for Elliptic Partial Differential Equations, 125196a0c994SBarry Smith Cambridge University Press. 1252b1a0a8a3SJed Brown 1253b1a0a8a3SJed Brown .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 125449517cdeSBarry Smith PCBJACOBI, PCGASMGetSubKSP(), PCGASMSetSubdomains(), 125534a84908Sprj- PCSetModifySubMatrices(), PCGASMSetOverlap(), PCGASMSetType() 1256b1a0a8a3SJed Brown 1257b1a0a8a3SJed Brown M*/ 1258b1a0a8a3SJed Brown 12598cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_GASM(PC pc) 1260b1a0a8a3SJed Brown { 1261b1a0a8a3SJed Brown PetscErrorCode ierr; 1262f746d493SDmitry Karpeev PC_GASM *osm; 1263b1a0a8a3SJed Brown 1264b1a0a8a3SJed Brown PetscFunctionBegin; 1265b00a9115SJed Brown ierr = PetscNewLog(pc,&osm);CHKERRQ(ierr); 12662fa5cd67SKarl Rupp 12678f3b4b4dSDmitry Karpeev osm->N = PETSC_DETERMINE; 126806b43e7eSDmitry Karpeev osm->n = PETSC_DECIDE; 12698f3b4b4dSDmitry Karpeev osm->nmax = PETSC_DETERMINE; 12708f3b4b4dSDmitry Karpeev osm->overlap = 0; 12710a545947SLisandro Dalcin osm->ksp = NULL; 12720a545947SLisandro Dalcin osm->gorestriction = NULL; 12730a545947SLisandro Dalcin osm->girestriction = NULL; 12740a545947SLisandro Dalcin osm->pctoouter = NULL; 12750a545947SLisandro Dalcin osm->gx = NULL; 12760a545947SLisandro Dalcin osm->gy = NULL; 12770a545947SLisandro Dalcin osm->x = NULL; 12780a545947SLisandro Dalcin osm->y = NULL; 12790a545947SLisandro Dalcin osm->pcx = NULL; 12800a545947SLisandro Dalcin osm->pcy = NULL; 12810a545947SLisandro Dalcin osm->permutationIS = NULL; 12820a545947SLisandro Dalcin osm->permutationP = NULL; 12830a545947SLisandro Dalcin osm->pcmat = NULL; 12840a545947SLisandro Dalcin osm->ois = NULL; 12850a545947SLisandro Dalcin osm->iis = NULL; 12860a545947SLisandro Dalcin osm->pmat = NULL; 1287f746d493SDmitry Karpeev osm->type = PC_GASM_RESTRICT; 12886a4f0f73SDmitry Karpeev osm->same_subdomain_solvers = PETSC_TRUE; 1289b1a0a8a3SJed Brown osm->sort_indices = PETSC_TRUE; 1290d709ea83SDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 1291ea91fabdSFande Kong osm->hierarchicalpartitioning = PETSC_FALSE; 1292b1a0a8a3SJed Brown 1293b1a0a8a3SJed Brown pc->data = (void*)osm; 1294f746d493SDmitry Karpeev pc->ops->apply = PCApply_GASM; 1295f746d493SDmitry Karpeev pc->ops->applytranspose = PCApplyTranspose_GASM; 1296f746d493SDmitry Karpeev pc->ops->setup = PCSetUp_GASM; 1297a06653b4SBarry Smith pc->ops->reset = PCReset_GASM; 1298f746d493SDmitry Karpeev pc->ops->destroy = PCDestroy_GASM; 1299f746d493SDmitry Karpeev pc->ops->setfromoptions = PCSetFromOptions_GASM; 1300f746d493SDmitry Karpeev pc->ops->setuponblocks = PCSetUpOnBlocks_GASM; 1301f746d493SDmitry Karpeev pc->ops->view = PCView_GASM; 13020a545947SLisandro Dalcin pc->ops->applyrichardson = NULL; 1303b1a0a8a3SJed Brown 1304bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetSubdomains_C",PCGASMSetSubdomains_GASM);CHKERRQ(ierr); 1305bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetOverlap_C",PCGASMSetOverlap_GASM);CHKERRQ(ierr); 1306bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetType_C",PCGASMSetType_GASM);CHKERRQ(ierr); 1307bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetSortIndices_C",PCGASMSetSortIndices_GASM);CHKERRQ(ierr); 1308bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGASMGetSubKSP_C",PCGASMGetSubKSP_GASM);CHKERRQ(ierr); 1309b1a0a8a3SJed Brown PetscFunctionReturn(0); 1310b1a0a8a3SJed Brown } 1311b1a0a8a3SJed Brown 1312b1a0a8a3SJed Brown 13138f3b4b4dSDmitry Karpeev PetscErrorCode PCGASMCreateLocalSubdomains(Mat A, PetscInt nloc, IS *iis[]) 1314b1a0a8a3SJed Brown { 1315b1a0a8a3SJed Brown MatPartitioning mpart; 1316b1a0a8a3SJed Brown const char *prefix; 1317b1a0a8a3SJed Brown PetscInt i,j,rstart,rend,bs; 1318976e8c5aSLisandro Dalcin PetscBool hasop, isbaij = PETSC_FALSE,foundpart = PETSC_FALSE; 13190298fd71SBarry Smith Mat Ad = NULL, adj; 1320b1a0a8a3SJed Brown IS ispart,isnumb,*is; 1321b1a0a8a3SJed Brown PetscErrorCode ierr; 1322b1a0a8a3SJed Brown 1323b1a0a8a3SJed Brown PetscFunctionBegin; 13248f3b4b4dSDmitry Karpeev if (nloc < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of local subdomains must > 0, got nloc = %D",nloc); 1325b1a0a8a3SJed Brown 1326b1a0a8a3SJed Brown /* Get prefix, row distribution, and block size */ 1327b1a0a8a3SJed Brown ierr = MatGetOptionsPrefix(A,&prefix);CHKERRQ(ierr); 1328b1a0a8a3SJed Brown ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 1329b1a0a8a3SJed Brown ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr); 1330b1a0a8a3SJed 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); 1331b1a0a8a3SJed Brown 1332b1a0a8a3SJed Brown /* Get diagonal block from matrix if possible */ 1333976e8c5aSLisandro Dalcin ierr = MatHasOperation(A,MATOP_GET_DIAGONAL_BLOCK,&hasop);CHKERRQ(ierr); 1334976e8c5aSLisandro Dalcin if (hasop) { 133511bd1e4dSLisandro Dalcin ierr = MatGetDiagonalBlock(A,&Ad);CHKERRQ(ierr); 1336b1a0a8a3SJed Brown } 1337b1a0a8a3SJed Brown if (Ad) { 1338b9e7e5c1SBarry Smith ierr = PetscObjectBaseTypeCompare((PetscObject)Ad,MATSEQBAIJ,&isbaij);CHKERRQ(ierr); 1339b9e7e5c1SBarry Smith if (!isbaij) {ierr = PetscObjectBaseTypeCompare((PetscObject)Ad,MATSEQSBAIJ,&isbaij);CHKERRQ(ierr);} 1340b1a0a8a3SJed Brown } 13418f3b4b4dSDmitry Karpeev if (Ad && nloc > 1) { 1342b1a0a8a3SJed Brown PetscBool match,done; 1343b1a0a8a3SJed Brown /* Try to setup a good matrix partitioning if available */ 1344b1a0a8a3SJed Brown ierr = MatPartitioningCreate(PETSC_COMM_SELF,&mpart);CHKERRQ(ierr); 1345b1a0a8a3SJed Brown ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr); 1346b1a0a8a3SJed Brown ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr); 1347251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGCURRENT,&match);CHKERRQ(ierr); 1348b1a0a8a3SJed Brown if (!match) { 1349251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGSQUARE,&match);CHKERRQ(ierr); 1350b1a0a8a3SJed Brown } 1351b1a0a8a3SJed Brown if (!match) { /* assume a "good" partitioner is available */ 13521a83f524SJed Brown PetscInt na; 13531a83f524SJed Brown const PetscInt *ia,*ja; 1354b1a0a8a3SJed Brown ierr = MatGetRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr); 1355b1a0a8a3SJed Brown if (done) { 1356b1a0a8a3SJed Brown /* Build adjacency matrix by hand. Unfortunately a call to 1357b1a0a8a3SJed Brown MatConvert(Ad,MATMPIADJ,MAT_INITIAL_MATRIX,&adj) will 1358b1a0a8a3SJed Brown remove the block-aij structure and we cannot expect 1359b1a0a8a3SJed Brown MatPartitioning to split vertices as we need */ 13600a545947SLisandro Dalcin PetscInt i,j,len,nnz,cnt,*iia=NULL,*jja=NULL; 13611a83f524SJed Brown const PetscInt *row; 1362b1a0a8a3SJed Brown nnz = 0; 1363b1a0a8a3SJed Brown for (i=0; i<na; i++) { /* count number of nonzeros */ 1364b1a0a8a3SJed Brown len = ia[i+1] - ia[i]; 1365b1a0a8a3SJed Brown row = ja + ia[i]; 1366b1a0a8a3SJed Brown for (j=0; j<len; j++) { 1367b1a0a8a3SJed Brown if (row[j] == i) { /* don't count diagonal */ 1368b1a0a8a3SJed Brown len--; break; 1369b1a0a8a3SJed Brown } 1370b1a0a8a3SJed Brown } 1371b1a0a8a3SJed Brown nnz += len; 1372b1a0a8a3SJed Brown } 1373854ce69bSBarry Smith ierr = PetscMalloc1(na+1,&iia);CHKERRQ(ierr); 1374854ce69bSBarry Smith ierr = PetscMalloc1(nnz,&jja);CHKERRQ(ierr); 1375b1a0a8a3SJed Brown nnz = 0; 1376b1a0a8a3SJed Brown iia[0] = 0; 1377b1a0a8a3SJed Brown for (i=0; i<na; i++) { /* fill adjacency */ 1378b1a0a8a3SJed Brown cnt = 0; 1379b1a0a8a3SJed Brown len = ia[i+1] - ia[i]; 1380b1a0a8a3SJed Brown row = ja + ia[i]; 1381b1a0a8a3SJed Brown for (j=0; j<len; j++) { 13822fa5cd67SKarl Rupp if (row[j] != i) jja[nnz+cnt++] = row[j]; /* if not diagonal */ 1383b1a0a8a3SJed Brown } 1384b1a0a8a3SJed Brown nnz += cnt; 1385b1a0a8a3SJed Brown iia[i+1] = nnz; 1386b1a0a8a3SJed Brown } 1387b1a0a8a3SJed Brown /* Partitioning of the adjacency matrix */ 13880298fd71SBarry Smith ierr = MatCreateMPIAdj(PETSC_COMM_SELF,na,na,iia,jja,NULL,&adj);CHKERRQ(ierr); 1389b1a0a8a3SJed Brown ierr = MatPartitioningSetAdjacency(mpart,adj);CHKERRQ(ierr); 13908f3b4b4dSDmitry Karpeev ierr = MatPartitioningSetNParts(mpart,nloc);CHKERRQ(ierr); 1391b1a0a8a3SJed Brown ierr = MatPartitioningApply(mpart,&ispart);CHKERRQ(ierr); 1392b1a0a8a3SJed Brown ierr = ISPartitioningToNumbering(ispart,&isnumb);CHKERRQ(ierr); 13936bf464f9SBarry Smith ierr = MatDestroy(&adj);CHKERRQ(ierr); 1394b1a0a8a3SJed Brown foundpart = PETSC_TRUE; 1395b1a0a8a3SJed Brown } 1396b1a0a8a3SJed Brown ierr = MatRestoreRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr); 1397b1a0a8a3SJed Brown } 13986bf464f9SBarry Smith ierr = MatPartitioningDestroy(&mpart);CHKERRQ(ierr); 1399b1a0a8a3SJed Brown } 14008f3b4b4dSDmitry Karpeev ierr = PetscMalloc1(nloc,&is);CHKERRQ(ierr); 1401b1a0a8a3SJed Brown if (!foundpart) { 1402b1a0a8a3SJed Brown 1403b1a0a8a3SJed Brown /* Partitioning by contiguous chunks of rows */ 1404b1a0a8a3SJed Brown 1405b1a0a8a3SJed Brown PetscInt mbs = (rend-rstart)/bs; 1406b1a0a8a3SJed Brown PetscInt start = rstart; 14078f3b4b4dSDmitry Karpeev for (i=0; i<nloc; i++) { 14088f3b4b4dSDmitry Karpeev PetscInt count = (mbs/nloc + ((mbs % nloc) > i)) * bs; 1409b1a0a8a3SJed Brown ierr = ISCreateStride(PETSC_COMM_SELF,count,start,1,&is[i]);CHKERRQ(ierr); 1410b1a0a8a3SJed Brown start += count; 1411b1a0a8a3SJed Brown } 1412b1a0a8a3SJed Brown 1413b1a0a8a3SJed Brown } else { 1414b1a0a8a3SJed Brown 1415b1a0a8a3SJed Brown /* Partitioning by adjacency of diagonal block */ 1416b1a0a8a3SJed Brown 1417b1a0a8a3SJed Brown const PetscInt *numbering; 1418b1a0a8a3SJed Brown PetscInt *count,nidx,*indices,*newidx,start=0; 1419b1a0a8a3SJed Brown /* Get node count in each partition */ 14208f3b4b4dSDmitry Karpeev ierr = PetscMalloc1(nloc,&count);CHKERRQ(ierr); 14218f3b4b4dSDmitry Karpeev ierr = ISPartitioningCount(ispart,nloc,count);CHKERRQ(ierr); 1422b1a0a8a3SJed Brown if (isbaij && bs > 1) { /* adjust for the block-aij case */ 14238f3b4b4dSDmitry Karpeev for (i=0; i<nloc; i++) count[i] *= bs; 1424b1a0a8a3SJed Brown } 1425b1a0a8a3SJed Brown /* Build indices from node numbering */ 1426b1a0a8a3SJed Brown ierr = ISGetLocalSize(isnumb,&nidx);CHKERRQ(ierr); 1427785e854fSJed Brown ierr = PetscMalloc1(nidx,&indices);CHKERRQ(ierr); 1428b1a0a8a3SJed Brown for (i=0; i<nidx; i++) indices[i] = i; /* needs to be initialized */ 1429b1a0a8a3SJed Brown ierr = ISGetIndices(isnumb,&numbering);CHKERRQ(ierr); 1430b1a0a8a3SJed Brown ierr = PetscSortIntWithPermutation(nidx,numbering,indices);CHKERRQ(ierr); 1431b1a0a8a3SJed Brown ierr = ISRestoreIndices(isnumb,&numbering);CHKERRQ(ierr); 1432b1a0a8a3SJed Brown if (isbaij && bs > 1) { /* adjust for the block-aij case */ 1433785e854fSJed Brown ierr = PetscMalloc1(nidx*bs,&newidx);CHKERRQ(ierr); 14342fa5cd67SKarl Rupp for (i=0; i<nidx; i++) { 14352fa5cd67SKarl Rupp for (j=0; j<bs; j++) newidx[i*bs+j] = indices[i]*bs + j; 14362fa5cd67SKarl Rupp } 1437b1a0a8a3SJed Brown ierr = PetscFree(indices);CHKERRQ(ierr); 1438b1a0a8a3SJed Brown nidx *= bs; 1439b1a0a8a3SJed Brown indices = newidx; 1440b1a0a8a3SJed Brown } 1441b1a0a8a3SJed Brown /* Shift to get global indices */ 1442b1a0a8a3SJed Brown for (i=0; i<nidx; i++) indices[i] += rstart; 1443b1a0a8a3SJed Brown 1444b1a0a8a3SJed Brown /* Build the index sets for each block */ 14458f3b4b4dSDmitry Karpeev for (i=0; i<nloc; i++) { 1446b1a0a8a3SJed Brown ierr = ISCreateGeneral(PETSC_COMM_SELF,count[i],&indices[start],PETSC_COPY_VALUES,&is[i]);CHKERRQ(ierr); 1447b1a0a8a3SJed Brown ierr = ISSort(is[i]);CHKERRQ(ierr); 1448b1a0a8a3SJed Brown start += count[i]; 1449b1a0a8a3SJed Brown } 1450b1a0a8a3SJed Brown 1451302440fdSBarry Smith ierr = PetscFree(count);CHKERRQ(ierr); 1452302440fdSBarry Smith ierr = PetscFree(indices);CHKERRQ(ierr); 1453fcfd50ebSBarry Smith ierr = ISDestroy(&isnumb);CHKERRQ(ierr); 1454fcfd50ebSBarry Smith ierr = ISDestroy(&ispart);CHKERRQ(ierr); 1455b1a0a8a3SJed Brown } 14566a4f0f73SDmitry Karpeev *iis = is; 14578f3b4b4dSDmitry Karpeev PetscFunctionReturn(0); 14588f3b4b4dSDmitry Karpeev } 14598f3b4b4dSDmitry Karpeev 1460b541e6a4SDmitry Karpeev PETSC_INTERN PetscErrorCode PCGASMCreateStraddlingSubdomains(Mat A,PetscInt N,PetscInt *n,IS *iis[]) 14618f3b4b4dSDmitry Karpeev { 14628f3b4b4dSDmitry Karpeev PetscErrorCode ierr; 14638f3b4b4dSDmitry Karpeev 14648f3b4b4dSDmitry Karpeev PetscFunctionBegin; 14658f3b4b4dSDmitry Karpeev ierr = MatSubdomainsCreateCoalesce(A,N,n,iis);CHKERRQ(ierr); 14668f3b4b4dSDmitry Karpeev PetscFunctionReturn(0); 14678f3b4b4dSDmitry Karpeev } 14688f3b4b4dSDmitry Karpeev 14698f3b4b4dSDmitry Karpeev 14708f3b4b4dSDmitry Karpeev 14718f3b4b4dSDmitry Karpeev /*@C 14728f3b4b4dSDmitry Karpeev PCGASMCreateSubdomains - Creates n index sets defining n nonoverlapping subdomains for the additive 14738f3b4b4dSDmitry Karpeev Schwarz preconditioner for a any problem based on its matrix. 14748f3b4b4dSDmitry Karpeev 14758f3b4b4dSDmitry Karpeev Collective 14768f3b4b4dSDmitry Karpeev 14778f3b4b4dSDmitry Karpeev Input Parameters: 14788f3b4b4dSDmitry Karpeev + A - The global matrix operator 14798f3b4b4dSDmitry Karpeev - N - the number of global subdomains requested 14808f3b4b4dSDmitry Karpeev 14818f3b4b4dSDmitry Karpeev Output Parameters: 14828f3b4b4dSDmitry Karpeev + n - the number of subdomains created on this processor 14838f3b4b4dSDmitry Karpeev - iis - the array of index sets defining the local inner subdomains (on which the correction is applied) 14848f3b4b4dSDmitry Karpeev 14858f3b4b4dSDmitry Karpeev Level: advanced 14868f3b4b4dSDmitry Karpeev 14878f3b4b4dSDmitry Karpeev Note: When N >= A's communicator size, each subdomain is local -- contained within a single processor. 14888f3b4b4dSDmitry Karpeev When N < size, the subdomains are 'straddling' (processor boundaries) and are no longer local. 14898f3b4b4dSDmitry Karpeev The resulting subdomains can be use in PCGASMSetSubdomains(pc,n,iss,NULL). The overlapping 14908f3b4b4dSDmitry Karpeev outer subdomains will be automatically generated from these according to the requested amount of 14918f3b4b4dSDmitry Karpeev overlap; this is currently supported only with local subdomains. 14928f3b4b4dSDmitry Karpeev 14938f3b4b4dSDmitry Karpeev 14948f3b4b4dSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMDestroySubdomains() 14958f3b4b4dSDmitry Karpeev @*/ 1496b541e6a4SDmitry Karpeev PetscErrorCode PCGASMCreateSubdomains(Mat A,PetscInt N,PetscInt *n,IS *iis[]) 14978f3b4b4dSDmitry Karpeev { 14988f3b4b4dSDmitry Karpeev PetscMPIInt size; 14998f3b4b4dSDmitry Karpeev PetscErrorCode ierr; 15008f3b4b4dSDmitry Karpeev 15018f3b4b4dSDmitry Karpeev PetscFunctionBegin; 15028f3b4b4dSDmitry Karpeev PetscValidHeaderSpecific(A,MAT_CLASSID,1); 15038f3b4b4dSDmitry Karpeev PetscValidPointer(iis,4); 15048f3b4b4dSDmitry Karpeev 15058f3b4b4dSDmitry Karpeev if (N < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Number of subdomains must be > 0, N = %D",N); 15068f3b4b4dSDmitry Karpeev ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr); 15078f3b4b4dSDmitry Karpeev if (N >= size) { 15088f3b4b4dSDmitry Karpeev *n = N/size + (N%size); 15098f3b4b4dSDmitry Karpeev ierr = PCGASMCreateLocalSubdomains(A,*n,iis);CHKERRQ(ierr); 15106a4f0f73SDmitry Karpeev } else { 15118f3b4b4dSDmitry Karpeev ierr = PCGASMCreateStraddlingSubdomains(A,N,n,iis);CHKERRQ(ierr); 15126a4f0f73SDmitry Karpeev } 1513b1a0a8a3SJed Brown PetscFunctionReturn(0); 1514b1a0a8a3SJed Brown } 1515b1a0a8a3SJed Brown 1516b1a0a8a3SJed Brown /*@C 1517f746d493SDmitry Karpeev PCGASMDestroySubdomains - Destroys the index sets created with 15188f3b4b4dSDmitry Karpeev PCGASMCreateSubdomains() or PCGASMCreateSubdomains2D. Should be 151906b43e7eSDmitry Karpeev called after setting subdomains with PCGASMSetSubdomains(). 1520b1a0a8a3SJed Brown 1521b1a0a8a3SJed Brown Collective 1522b1a0a8a3SJed Brown 1523b1a0a8a3SJed Brown Input Parameters: 1524b1a0a8a3SJed Brown + n - the number of index sets 15256a4f0f73SDmitry Karpeev . iis - the array of inner subdomains, 15260298fd71SBarry Smith - ois - the array of outer subdomains, can be NULL 1527b1a0a8a3SJed Brown 15286a4f0f73SDmitry Karpeev Level: intermediate 15296a4f0f73SDmitry Karpeev 153095452b02SPatrick Sanan Notes: 153195452b02SPatrick Sanan this is merely a convenience subroutine that walks each list, 15322c112581SDmitry Karpeev destroys each IS on the list, and then frees the list. At the end the 15332c112581SDmitry Karpeev list pointers are set to NULL. 1534b1a0a8a3SJed Brown 15358f3b4b4dSDmitry Karpeev .seealso: PCGASMCreateSubdomains(), PCGASMSetSubdomains() 1536b1a0a8a3SJed Brown @*/ 15372c112581SDmitry Karpeev PetscErrorCode PCGASMDestroySubdomains(PetscInt n,IS **iis,IS **ois) 1538b1a0a8a3SJed Brown { 1539b1a0a8a3SJed Brown PetscInt i; 1540b1a0a8a3SJed Brown PetscErrorCode ierr; 15415fd66863SKarl Rupp 1542b1a0a8a3SJed Brown PetscFunctionBegin; 1543a35b7b57SDmitry Karpeev if (n <= 0) PetscFunctionReturn(0); 15446a4f0f73SDmitry Karpeev if (ois) { 15452c112581SDmitry Karpeev PetscValidPointer(ois,3); 15462c112581SDmitry Karpeev if (*ois) { 15472c112581SDmitry Karpeev PetscValidPointer(*ois,3); 1548a35b7b57SDmitry Karpeev for (i=0; i<n; i++) { 15492c112581SDmitry Karpeev ierr = ISDestroy(&(*ois)[i]);CHKERRQ(ierr); 1550a35b7b57SDmitry Karpeev } 1551135757f6SDmitry Karpeev ierr = PetscFree((*ois));CHKERRQ(ierr); 15522c112581SDmitry Karpeev } 1553b1a0a8a3SJed Brown } 1554b9d0fdaaSFande Kong if (iis) { 1555b9d0fdaaSFande Kong PetscValidPointer(iis,2); 1556b9d0fdaaSFande Kong if (*iis) { 1557b9d0fdaaSFande Kong PetscValidPointer(*iis,2); 1558b9d0fdaaSFande Kong for (i=0; i<n; i++) { 1559b9d0fdaaSFande Kong ierr = ISDestroy(&(*iis)[i]);CHKERRQ(ierr); 1560b9d0fdaaSFande Kong } 1561b9d0fdaaSFande Kong ierr = PetscFree((*iis));CHKERRQ(ierr); 1562b9d0fdaaSFande Kong } 1563b9d0fdaaSFande Kong } 1564b1a0a8a3SJed Brown PetscFunctionReturn(0); 1565b1a0a8a3SJed Brown } 1566b1a0a8a3SJed Brown 1567af538404SDmitry Karpeev 1568af538404SDmitry Karpeev #define PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,xleft_loc,ylow_loc,xright_loc,yhigh_loc,n) \ 1569af538404SDmitry Karpeev { \ 1570af538404SDmitry Karpeev PetscInt first_row = first/M, last_row = last/M+1; \ 1571af538404SDmitry Karpeev /* \ 1572af538404SDmitry Karpeev Compute ylow_loc and yhigh_loc so that (ylow_loc,xleft) and (yhigh_loc,xright) are the corners \ 1573af538404SDmitry Karpeev of the bounding box of the intersection of the subdomain with the local ownership range (local \ 1574af538404SDmitry Karpeev subdomain). \ 1575af538404SDmitry Karpeev Also compute xleft_loc and xright_loc as the lower and upper bounds on the first and last rows \ 1576af538404SDmitry Karpeev of the intersection. \ 1577af538404SDmitry Karpeev */ \ 1578af538404SDmitry Karpeev /* ylow_loc is the grid row containing the first element of the local sumbdomain */ \ 1579eec7e3faSDmitry Karpeev *ylow_loc = PetscMax(first_row,ylow); \ 1580af538404SDmitry Karpeev /* xleft_loc is the offset of first element of the local subdomain within its grid row (might actually be outside the local subdomain) */ \ 1581eec7e3faSDmitry Karpeev *xleft_loc = *ylow_loc==first_row ? PetscMax(first%M,xleft) : xleft; \ 1582af538404SDmitry Karpeev /* yhigh_loc is the grid row above the last local subdomain element */ \ 1583eec7e3faSDmitry Karpeev *yhigh_loc = PetscMin(last_row,yhigh); \ 1584af538404SDmitry Karpeev /* xright is the offset of the end of the local subdomain within its grid row (might actually be outside the local subdomain) */ \ 1585eec7e3faSDmitry Karpeev *xright_loc = *yhigh_loc==last_row ? PetscMin(xright,last%M) : xright; \ 1586af538404SDmitry Karpeev /* Now compute the size of the local subdomain n. */ \ 1587c3518bceSDmitry Karpeev *n = 0; \ 1588eec7e3faSDmitry Karpeev if (*ylow_loc < *yhigh_loc) { \ 1589af538404SDmitry Karpeev PetscInt width = xright-xleft; \ 1590eec7e3faSDmitry Karpeev *n += width*(*yhigh_loc-*ylow_loc-1); \ 1591eec7e3faSDmitry Karpeev *n += PetscMin(PetscMax(*xright_loc-xleft,0),width); \ 1592eec7e3faSDmitry Karpeev *n -= PetscMin(PetscMax(*xleft_loc-xleft,0), width); \ 1593af538404SDmitry Karpeev } \ 1594af538404SDmitry Karpeev } 1595af538404SDmitry Karpeev 1596c3518bceSDmitry Karpeev 1597eec7e3faSDmitry Karpeev 1598b1a0a8a3SJed Brown /*@ 1599f746d493SDmitry Karpeev PCGASMCreateSubdomains2D - Creates the index sets for the overlapping Schwarz 1600b1a0a8a3SJed Brown preconditioner for a two-dimensional problem on a regular grid. 1601b1a0a8a3SJed Brown 1602af538404SDmitry Karpeev Collective 1603b1a0a8a3SJed Brown 1604b1a0a8a3SJed Brown Input Parameters: 160506b43e7eSDmitry Karpeev + M, N - the global number of grid points in the x and y directions 1606af538404SDmitry Karpeev . Mdomains, Ndomains - the global number of subdomains in the x and y directions 1607b1a0a8a3SJed Brown . dof - degrees of freedom per node 1608b1a0a8a3SJed Brown - overlap - overlap in mesh lines 1609b1a0a8a3SJed Brown 1610b1a0a8a3SJed Brown Output Parameters: 1611af538404SDmitry Karpeev + Nsub - the number of local subdomains created 16126a4f0f73SDmitry Karpeev . iis - array of index sets defining inner (nonoverlapping) subdomains 16136a4f0f73SDmitry Karpeev - ois - array of index sets defining outer (overlapping, if overlap > 0) subdomains 1614b1a0a8a3SJed Brown 1615b1a0a8a3SJed Brown 1616b1a0a8a3SJed Brown Level: advanced 1617b1a0a8a3SJed Brown 16188f3b4b4dSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMGetSubKSP(), PCGASMSetOverlap() 1619b1a0a8a3SJed Brown @*/ 16206a4f0f73SDmitry Karpeev PetscErrorCode PCGASMCreateSubdomains2D(PC pc,PetscInt M,PetscInt N,PetscInt Mdomains,PetscInt Ndomains,PetscInt dof,PetscInt overlap,PetscInt *nsub,IS **iis,IS **ois) 1621b1a0a8a3SJed Brown { 1622b1a0a8a3SJed Brown PetscErrorCode ierr; 16232bbb417fSDmitry Karpeev PetscMPIInt size, rank; 16242bbb417fSDmitry Karpeev PetscInt i, j; 16252bbb417fSDmitry Karpeev PetscInt maxheight, maxwidth; 16262bbb417fSDmitry Karpeev PetscInt xstart, xleft, xright, xleft_loc, xright_loc; 16272bbb417fSDmitry Karpeev PetscInt ystart, ylow, yhigh, ylow_loc, yhigh_loc; 1628eec7e3faSDmitry Karpeev PetscInt x[2][2], y[2][2], n[2]; 16292bbb417fSDmitry Karpeev PetscInt first, last; 16302bbb417fSDmitry Karpeev PetscInt nidx, *idx; 16312bbb417fSDmitry Karpeev PetscInt ii,jj,s,q,d; 1632af538404SDmitry Karpeev PetscInt k,kk; 16332bbb417fSDmitry Karpeev PetscMPIInt color; 16342bbb417fSDmitry Karpeev MPI_Comm comm, subcomm; 16350a545947SLisandro Dalcin IS **xis = NULL, **is = ois, **is_local = iis; 1636d34fcf5fSBarry Smith 1637b1a0a8a3SJed Brown PetscFunctionBegin; 16382bbb417fSDmitry Karpeev ierr = PetscObjectGetComm((PetscObject)pc, &comm);CHKERRQ(ierr); 16392bbb417fSDmitry Karpeev ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 16402bbb417fSDmitry Karpeev ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 16412bbb417fSDmitry Karpeev ierr = MatGetOwnershipRange(pc->pmat, &first, &last);CHKERRQ(ierr); 1642d34fcf5fSBarry 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) " 16432bbb417fSDmitry Karpeev "does not respect the number of degrees of freedom per grid point %D", first, last, dof); 1644d34fcf5fSBarry Smith 1645af538404SDmitry Karpeev /* Determine the number of domains with nonzero intersections with the local ownership range. */ 16462bbb417fSDmitry Karpeev s = 0; 1647af538404SDmitry Karpeev ystart = 0; 1648af538404SDmitry Karpeev for (j=0; j<Ndomains; ++j) { 1649af538404SDmitry Karpeev maxheight = N/Ndomains + ((N % Ndomains) > j); /* Maximal height of subdomain */ 1650af538404SDmitry 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); 1651eec7e3faSDmitry Karpeev /* Vertical domain limits with an overlap. */ 1652eec7e3faSDmitry Karpeev ylow = PetscMax(ystart - overlap,0); 1653eec7e3faSDmitry Karpeev yhigh = PetscMin(ystart + maxheight + overlap,N); 1654b1a0a8a3SJed Brown xstart = 0; 1655af538404SDmitry Karpeev for (i=0; i<Mdomains; ++i) { 1656af538404SDmitry Karpeev maxwidth = M/Mdomains + ((M % Mdomains) > i); /* Maximal width of subdomain */ 1657af538404SDmitry 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); 1658eec7e3faSDmitry Karpeev /* Horizontal domain limits with an overlap. */ 1659eec7e3faSDmitry Karpeev xleft = PetscMax(xstart - overlap,0); 1660eec7e3faSDmitry Karpeev xright = PetscMin(xstart + maxwidth + overlap,M); 1661af538404SDmitry Karpeev /* 1662af538404SDmitry Karpeev Determine whether this subdomain intersects this processor's ownership range of pc->pmat. 1663af538404SDmitry Karpeev */ 1664c3518bceSDmitry Karpeev PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,(&xleft_loc),(&ylow_loc),(&xright_loc),(&yhigh_loc),(&nidx)); 16652fa5cd67SKarl Rupp if (nidx) ++s; 1666af538404SDmitry Karpeev xstart += maxwidth; 1667af538404SDmitry Karpeev } /* for (i = 0; i < Mdomains; ++i) */ 1668af538404SDmitry Karpeev ystart += maxheight; 1669af538404SDmitry Karpeev } /* for (j = 0; j < Ndomains; ++j) */ 16702fa5cd67SKarl Rupp 1671af538404SDmitry Karpeev /* Now we can allocate the necessary number of ISs. */ 1672af538404SDmitry Karpeev *nsub = s; 1673854ce69bSBarry Smith ierr = PetscMalloc1(*nsub,is);CHKERRQ(ierr); 1674854ce69bSBarry Smith ierr = PetscMalloc1(*nsub,is_local);CHKERRQ(ierr); 1675af538404SDmitry Karpeev s = 0; 1676af538404SDmitry Karpeev ystart = 0; 1677af538404SDmitry Karpeev for (j=0; j<Ndomains; ++j) { 1678af538404SDmitry Karpeev maxheight = N/Ndomains + ((N % Ndomains) > j); /* Maximal height of subdomain */ 1679af538404SDmitry 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); 1680af538404SDmitry Karpeev /* Vertical domain limits with an overlap. */ 1681eec7e3faSDmitry Karpeev y[0][0] = PetscMax(ystart - overlap,0); 1682eec7e3faSDmitry Karpeev y[0][1] = PetscMin(ystart + maxheight + overlap,N); 1683eec7e3faSDmitry Karpeev /* Vertical domain limits without an overlap. */ 1684eec7e3faSDmitry Karpeev y[1][0] = ystart; 1685eec7e3faSDmitry Karpeev y[1][1] = PetscMin(ystart + maxheight,N); 1686eec7e3faSDmitry Karpeev xstart = 0; 1687af538404SDmitry Karpeev for (i=0; i<Mdomains; ++i) { 1688af538404SDmitry Karpeev maxwidth = M/Mdomains + ((M % Mdomains) > i); /* Maximal width of subdomain */ 1689af538404SDmitry 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); 1690af538404SDmitry Karpeev /* Horizontal domain limits with an overlap. */ 1691eec7e3faSDmitry Karpeev x[0][0] = PetscMax(xstart - overlap,0); 1692eec7e3faSDmitry Karpeev x[0][1] = PetscMin(xstart + maxwidth + overlap,M); 1693eec7e3faSDmitry Karpeev /* Horizontal domain limits without an overlap. */ 1694eec7e3faSDmitry Karpeev x[1][0] = xstart; 1695eec7e3faSDmitry Karpeev x[1][1] = PetscMin(xstart+maxwidth,M); 16962bbb417fSDmitry Karpeev /* 16972bbb417fSDmitry Karpeev Determine whether this domain intersects this processor's ownership range of pc->pmat. 16982bbb417fSDmitry Karpeev Do this twice: first for the domains with overlaps, and once without. 16992bbb417fSDmitry Karpeev During the first pass create the subcommunicators, and use them on the second pass as well. 17002bbb417fSDmitry Karpeev */ 17012bbb417fSDmitry Karpeev for (q = 0; q < 2; ++q) { 1702cc96b2e5SBarry Smith PetscBool split = PETSC_FALSE; 17032bbb417fSDmitry Karpeev /* 17042bbb417fSDmitry Karpeev domain limits, (xleft, xright) and (ylow, yheigh) are adjusted 17052bbb417fSDmitry Karpeev according to whether the domain with an overlap or without is considered. 17062bbb417fSDmitry Karpeev */ 17072bbb417fSDmitry Karpeev xleft = x[q][0]; xright = x[q][1]; 17082bbb417fSDmitry Karpeev ylow = y[q][0]; yhigh = y[q][1]; 1709c3518bceSDmitry Karpeev PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,(&xleft_loc),(&ylow_loc),(&xright_loc),(&yhigh_loc),(&nidx)); 1710af538404SDmitry Karpeev nidx *= dof; 1711eec7e3faSDmitry Karpeev n[q] = nidx; 17122bbb417fSDmitry Karpeev /* 1713eec7e3faSDmitry Karpeev Based on the counted number of indices in the local domain *with an overlap*, 1714af538404SDmitry Karpeev construct a subcommunicator of all the processors supporting this domain. 1715eec7e3faSDmitry Karpeev Observe that a domain with an overlap might have nontrivial local support, 1716eec7e3faSDmitry Karpeev while the domain without an overlap might not. Hence, the decision to participate 1717eec7e3faSDmitry Karpeev in the subcommunicator must be based on the domain with an overlap. 17182bbb417fSDmitry Karpeev */ 17192bbb417fSDmitry Karpeev if (q == 0) { 17202fa5cd67SKarl Rupp if (nidx) color = 1; 17212fa5cd67SKarl Rupp else color = MPI_UNDEFINED; 17222bbb417fSDmitry Karpeev ierr = MPI_Comm_split(comm, color, rank, &subcomm);CHKERRQ(ierr); 1723cc96b2e5SBarry Smith split = PETSC_TRUE; 17242bbb417fSDmitry Karpeev } 1725af538404SDmitry Karpeev /* 1726eec7e3faSDmitry Karpeev Proceed only if the number of local indices *with an overlap* is nonzero. 1727af538404SDmitry Karpeev */ 1728eec7e3faSDmitry Karpeev if (n[0]) { 17292fa5cd67SKarl Rupp if (q == 0) xis = is; 1730af538404SDmitry Karpeev if (q == 1) { 1731af538404SDmitry Karpeev /* 1732eec7e3faSDmitry Karpeev The IS for the no-overlap subdomain shares a communicator with the overlapping domain. 1733af538404SDmitry Karpeev Moreover, if the overlap is zero, the two ISs are identical. 1734af538404SDmitry Karpeev */ 1735b1a0a8a3SJed Brown if (overlap == 0) { 1736eec7e3faSDmitry Karpeev (*is_local)[s] = (*is)[s]; 17372bbb417fSDmitry Karpeev ierr = PetscObjectReference((PetscObject)(*is)[s]);CHKERRQ(ierr); 17382bbb417fSDmitry Karpeev continue; 1739d34fcf5fSBarry Smith } else { 17406a4f0f73SDmitry Karpeev xis = is_local; 1741eec7e3faSDmitry Karpeev subcomm = ((PetscObject)(*is)[s])->comm; 17422bbb417fSDmitry Karpeev } 1743af538404SDmitry Karpeev } /* if (q == 1) */ 17440298fd71SBarry Smith idx = NULL; 1745785e854fSJed Brown ierr = PetscMalloc1(nidx,&idx);CHKERRQ(ierr); 1746eec7e3faSDmitry Karpeev if (nidx) { 17472bbb417fSDmitry Karpeev k = 0; 17482bbb417fSDmitry Karpeev for (jj=ylow_loc; jj<yhigh_loc; ++jj) { 1749af538404SDmitry Karpeev PetscInt x0 = (jj==ylow_loc) ? xleft_loc : xleft; 1750af538404SDmitry Karpeev PetscInt x1 = (jj==yhigh_loc-1) ? xright_loc : xright; 1751af538404SDmitry Karpeev kk = dof*(M*jj + x0); 17522bbb417fSDmitry Karpeev for (ii=x0; ii<x1; ++ii) { 17532bbb417fSDmitry Karpeev for (d = 0; d < dof; ++d) { 17542bbb417fSDmitry Karpeev idx[k++] = kk++; 1755b1a0a8a3SJed Brown } 1756b1a0a8a3SJed Brown } 1757b1a0a8a3SJed Brown } 1758eec7e3faSDmitry Karpeev } 17596a4f0f73SDmitry Karpeev ierr = ISCreateGeneral(subcomm,nidx,idx,PETSC_OWN_POINTER,(*xis)+s);CHKERRQ(ierr); 1760cc96b2e5SBarry Smith if (split) { 1761cc96b2e5SBarry Smith ierr = MPI_Comm_free(&subcomm);CHKERRQ(ierr); 1762cc96b2e5SBarry Smith } 1763eec7e3faSDmitry Karpeev }/* if (n[0]) */ 17642bbb417fSDmitry Karpeev }/* for (q = 0; q < 2; ++q) */ 17652fa5cd67SKarl Rupp if (n[0]) ++s; 1766af538404SDmitry Karpeev xstart += maxwidth; 1767af538404SDmitry Karpeev } /* for (i = 0; i < Mdomains; ++i) */ 17682bbb417fSDmitry Karpeev ystart += maxheight; 1769af538404SDmitry Karpeev } /* for (j = 0; j < Ndomains; ++j) */ 1770b1a0a8a3SJed Brown PetscFunctionReturn(0); 1771b1a0a8a3SJed Brown } 1772b1a0a8a3SJed Brown 1773b1a0a8a3SJed Brown /*@C 177406b43e7eSDmitry Karpeev PCGASMGetSubdomains - Gets the subdomains supported on this processor 177506b43e7eSDmitry Karpeev for the additive Schwarz preconditioner. 1776b1a0a8a3SJed Brown 1777b1a0a8a3SJed Brown Not Collective 1778b1a0a8a3SJed Brown 1779b1a0a8a3SJed Brown Input Parameter: 1780b1a0a8a3SJed Brown . pc - the preconditioner context 1781b1a0a8a3SJed Brown 1782b1a0a8a3SJed Brown Output Parameters: 1783b1a0a8a3SJed Brown + n - the number of subdomains for this processor (default value = 1) 17840298fd71SBarry Smith . iis - the index sets that define the inner subdomains (without overlap) supported on this processor (can be NULL) 17850298fd71SBarry Smith - ois - the index sets that define the outer subdomains (with overlap) supported on this processor (can be NULL) 1786b1a0a8a3SJed Brown 1787b1a0a8a3SJed Brown 1788b1a0a8a3SJed Brown Notes: 17896a4f0f73SDmitry Karpeev The user is responsible for destroying the ISs and freeing the returned arrays. 1790b1a0a8a3SJed Brown The IS numbering is in the parallel, global numbering of the vector. 1791b1a0a8a3SJed Brown 1792b1a0a8a3SJed Brown Level: advanced 1793b1a0a8a3SJed Brown 17948f3b4b4dSDmitry Karpeev .seealso: PCGASMSetOverlap(), PCGASMGetSubKSP(), PCGASMCreateSubdomains2D(), 17958f3b4b4dSDmitry Karpeev PCGASMSetSubdomains(), PCGASMGetSubmatrices() 1796b1a0a8a3SJed Brown @*/ 17976a4f0f73SDmitry Karpeev PetscErrorCode PCGASMGetSubdomains(PC pc,PetscInt *n,IS *iis[],IS *ois[]) 1798b1a0a8a3SJed Brown { 1799f746d493SDmitry Karpeev PC_GASM *osm; 1800b1a0a8a3SJed Brown PetscErrorCode ierr; 1801b1a0a8a3SJed Brown PetscBool match; 18026a4f0f73SDmitry Karpeev PetscInt i; 18035fd66863SKarl Rupp 1804b1a0a8a3SJed Brown PetscFunctionBegin; 1805b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1806251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match);CHKERRQ(ierr); 1807ce94432eSBarry Smith if (!match) SETERRQ2(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Incorrect object type: expected %s, got %s instead", PCGASM, ((PetscObject)pc)->type_name); 1808f746d493SDmitry Karpeev osm = (PC_GASM*)pc->data; 1809ab3e923fSDmitry Karpeev if (n) *n = osm->n; 18106a4f0f73SDmitry Karpeev if (iis) { 1811785e854fSJed Brown ierr = PetscMalloc1(osm->n, iis);CHKERRQ(ierr); 18126a4f0f73SDmitry Karpeev } 18136a4f0f73SDmitry Karpeev if (ois) { 1814785e854fSJed Brown ierr = PetscMalloc1(osm->n, ois);CHKERRQ(ierr); 18156a4f0f73SDmitry Karpeev } 18166a4f0f73SDmitry Karpeev if (iis || ois) { 18176a4f0f73SDmitry Karpeev for (i = 0; i < osm->n; ++i) { 18186a4f0f73SDmitry Karpeev if (iis) (*iis)[i] = osm->iis[i]; 18196a4f0f73SDmitry Karpeev if (ois) (*ois)[i] = osm->ois[i]; 18206a4f0f73SDmitry Karpeev } 1821b1a0a8a3SJed Brown } 1822b1a0a8a3SJed Brown PetscFunctionReturn(0); 1823b1a0a8a3SJed Brown } 1824b1a0a8a3SJed Brown 1825b1a0a8a3SJed Brown /*@C 182606b43e7eSDmitry Karpeev PCGASMGetSubmatrices - Gets the local submatrices (for this processor 1827b1a0a8a3SJed Brown only) for the additive Schwarz preconditioner. 1828b1a0a8a3SJed Brown 1829b1a0a8a3SJed Brown Not Collective 1830b1a0a8a3SJed Brown 1831b1a0a8a3SJed Brown Input Parameter: 1832b1a0a8a3SJed Brown . pc - the preconditioner context 1833b1a0a8a3SJed Brown 1834b1a0a8a3SJed Brown Output Parameters: 1835b1a0a8a3SJed Brown + n - the number of matrices for this processor (default value = 1) 1836b1a0a8a3SJed Brown - mat - the matrices 1837b1a0a8a3SJed Brown 183895452b02SPatrick Sanan Notes: 183995452b02SPatrick Sanan matrices returned by this routine have the same communicators as the index sets (IS) 18408f3b4b4dSDmitry Karpeev used to define subdomains in PCGASMSetSubdomains() 1841b1a0a8a3SJed Brown Level: advanced 1842b1a0a8a3SJed Brown 18438f3b4b4dSDmitry Karpeev .seealso: PCGASMSetOverlap(), PCGASMGetSubKSP(), 184406b43e7eSDmitry Karpeev PCGASMCreateSubdomains2D(), PCGASMSetSubdomains(), PCGASMGetSubdomains() 1845b1a0a8a3SJed Brown @*/ 184606b43e7eSDmitry Karpeev PetscErrorCode PCGASMGetSubmatrices(PC pc,PetscInt *n,Mat *mat[]) 1847b1a0a8a3SJed Brown { 1848f746d493SDmitry Karpeev PC_GASM *osm; 1849b1a0a8a3SJed Brown PetscErrorCode ierr; 1850b1a0a8a3SJed Brown PetscBool match; 1851b1a0a8a3SJed Brown 1852b1a0a8a3SJed Brown PetscFunctionBegin; 1853b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1854b1a0a8a3SJed Brown PetscValidIntPointer(n,2); 1855b1a0a8a3SJed Brown if (mat) PetscValidPointer(mat,3); 185634a84908Sprj- if (!pc->setupcalled) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call after KSPSetUp() or PCSetUp()."); 1857251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match);CHKERRQ(ierr); 1858ce94432eSBarry Smith if (!match) SETERRQ2(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Expected %s, got %s instead", PCGASM, ((PetscObject)pc)->type_name); 1859f746d493SDmitry Karpeev osm = (PC_GASM*)pc->data; 1860ab3e923fSDmitry Karpeev if (n) *n = osm->n; 1861b1a0a8a3SJed Brown if (mat) *mat = osm->pmat; 1862b1a0a8a3SJed Brown PetscFunctionReturn(0); 1863b1a0a8a3SJed Brown } 1864d709ea83SDmitry Karpeev 1865d709ea83SDmitry Karpeev /*@ 18668f3b4b4dSDmitry Karpeev PCGASMSetUseDMSubdomains - Indicates whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible. 1867d709ea83SDmitry Karpeev Logically Collective 1868d709ea83SDmitry Karpeev 1869d709ea83SDmitry Karpeev Input Parameter: 1870d709ea83SDmitry Karpeev + pc - the preconditioner 1871d709ea83SDmitry Karpeev - flg - boolean indicating whether to use subdomains defined by the DM 1872d709ea83SDmitry Karpeev 1873d709ea83SDmitry Karpeev Options Database Key: 18748f3b4b4dSDmitry Karpeev . -pc_gasm_dm_subdomains -pc_gasm_overlap -pc_gasm_total_subdomains 1875d709ea83SDmitry Karpeev 1876d709ea83SDmitry Karpeev Level: intermediate 1877d709ea83SDmitry Karpeev 1878d709ea83SDmitry Karpeev Notes: 18798f3b4b4dSDmitry Karpeev PCGASMSetSubdomains(), PCGASMSetTotalSubdomains() or PCGASMSetOverlap() take precedence over PCGASMSetUseDMSubdomains(), 18808f3b4b4dSDmitry Karpeev so setting PCGASMSetSubdomains() with nontrivial subdomain ISs or any of PCGASMSetTotalSubdomains() and PCGASMSetOverlap() 18818f3b4b4dSDmitry Karpeev automatically turns the latter off. 1882d709ea83SDmitry Karpeev 18838f3b4b4dSDmitry Karpeev .seealso: PCGASMGetUseDMSubdomains(), PCGASMSetSubdomains(), PCGASMSetOverlap() 1884d709ea83SDmitry Karpeev PCGASMCreateSubdomains2D() 1885d709ea83SDmitry Karpeev @*/ 18868f3b4b4dSDmitry Karpeev PetscErrorCode PCGASMSetUseDMSubdomains(PC pc,PetscBool flg) 1887d709ea83SDmitry Karpeev { 1888d709ea83SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 1889d709ea83SDmitry Karpeev PetscErrorCode ierr; 1890d709ea83SDmitry Karpeev PetscBool match; 1891d709ea83SDmitry Karpeev 1892d709ea83SDmitry Karpeev PetscFunctionBegin; 1893d709ea83SDmitry Karpeev PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1894d709ea83SDmitry Karpeev PetscValidLogicalCollectiveBool(pc,flg,2); 1895d709ea83SDmitry Karpeev if (pc->setupcalled) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Not for a setup PC."); 1896d709ea83SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match);CHKERRQ(ierr); 1897d709ea83SDmitry Karpeev if (match) { 18988f3b4b4dSDmitry Karpeev if (!osm->user_subdomains && osm->N == PETSC_DETERMINE && osm->overlap < 0) { 1899d709ea83SDmitry Karpeev osm->dm_subdomains = flg; 1900d709ea83SDmitry Karpeev } 19018f3b4b4dSDmitry Karpeev } 1902d709ea83SDmitry Karpeev PetscFunctionReturn(0); 1903d709ea83SDmitry Karpeev } 1904d709ea83SDmitry Karpeev 1905d709ea83SDmitry Karpeev /*@ 19068f3b4b4dSDmitry Karpeev PCGASMGetUseDMSubdomains - Returns flag indicating whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible. 1907d709ea83SDmitry Karpeev Not Collective 1908d709ea83SDmitry Karpeev 1909d709ea83SDmitry Karpeev Input Parameter: 1910d709ea83SDmitry Karpeev . pc - the preconditioner 1911d709ea83SDmitry Karpeev 1912d709ea83SDmitry Karpeev Output Parameter: 1913d709ea83SDmitry Karpeev . flg - boolean indicating whether to use subdomains defined by the DM 1914d709ea83SDmitry Karpeev 1915d709ea83SDmitry Karpeev Level: intermediate 1916d709ea83SDmitry Karpeev 19178f3b4b4dSDmitry Karpeev .seealso: PCGASMSetUseDMSubdomains(), PCGASMSetOverlap() 1918d709ea83SDmitry Karpeev PCGASMCreateSubdomains2D() 1919d709ea83SDmitry Karpeev @*/ 19208f3b4b4dSDmitry Karpeev PetscErrorCode PCGASMGetUseDMSubdomains(PC pc,PetscBool* flg) 1921d709ea83SDmitry Karpeev { 1922d709ea83SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 1923d709ea83SDmitry Karpeev PetscErrorCode ierr; 1924d709ea83SDmitry Karpeev PetscBool match; 1925d709ea83SDmitry Karpeev 1926d709ea83SDmitry Karpeev PetscFunctionBegin; 1927d709ea83SDmitry Karpeev PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1928534a8f05SLisandro Dalcin PetscValidBoolPointer(flg,2); 1929d709ea83SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match);CHKERRQ(ierr); 1930d709ea83SDmitry Karpeev if (match) { 1931d709ea83SDmitry Karpeev if (flg) *flg = osm->dm_subdomains; 1932d709ea83SDmitry Karpeev } 1933d709ea83SDmitry Karpeev PetscFunctionReturn(0); 1934d709ea83SDmitry Karpeev } 1935