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 */ 73*36a9e3b9SBarry Smith #define len 16*(nidx+1)+1 74*36a9e3b9SBarry Smith ierr = PetscMalloc1(len, &cidx);CHKERRQ(ierr); 75*36a9e3b9SBarry Smith ierr = PetscViewerStringOpen(PETSC_COMM_SELF, cidx, len, &sviewer);CHKERRQ(ierr); 76*36a9e3b9SBarry 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 */ 100*36a9e3b9SBarry Smith #define len 16*(nidx+1)+1 101*36a9e3b9SBarry Smith ierr = PetscMalloc1(len, &cidx);CHKERRQ(ierr); 102*36a9e3b9SBarry Smith ierr = PetscViewerStringOpen(PETSC_COMM_SELF, cidx, len, &sviewer);CHKERRQ(ierr); 103*36a9e3b9SBarry 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); 138c5929fdfSBarry Smith ierr = PetscOptionsGetString(NULL,prefix,"-pc_gasm_print_subdomains",fname,PETSC_MAX_PATH_LEN,&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; 322b1a0a8a3SJed Brown PetscBool symset,flg; 323ea91fabdSFande Kong PetscInt i,nInnerIndices,nTotalInnerIndices; 324c10bc1a1SDmitry Karpeev PetscMPIInt rank, size; 325b1a0a8a3SJed Brown MatReuse scall = MAT_REUSE_MATRIX; 326b1a0a8a3SJed Brown KSP ksp; 327b1a0a8a3SJed Brown PC subpc; 328b1a0a8a3SJed Brown const char *prefix,*pprefix; 329f746d493SDmitry Karpeev Vec x,y; 3306a4f0f73SDmitry Karpeev PetscInt oni; /* Number of indices in the i-th local outer subdomain. */ 3316a4f0f73SDmitry Karpeev const PetscInt *oidxi; /* Indices from the i-th subdomain local outer subdomain. */ 3326a4f0f73SDmitry Karpeev PetscInt on; /* Number of indices in the disjoint union of local outer subdomains. */ 3336a4f0f73SDmitry Karpeev PetscInt *oidx; /* Indices in the disjoint union of local outer subdomains. */ 3346a4f0f73SDmitry Karpeev IS gois; /* Disjoint union the global indices of outer subdomains. */ 3356a4f0f73SDmitry Karpeev IS goid; /* Identity IS of the size of the disjoint union of outer subdomains. */ 336f746d493SDmitry Karpeev PetscScalar *gxarray, *gyarray; 337930d09c1SFande Kong PetscInt gostart; /* Start of locally-owned indices in the vectors -- osm->gx,osm->gy -- over the disjoint union of outer subdomains. */ 3388f3b4b4dSDmitry Karpeev PetscInt num_subdomains = 0; 3390298fd71SBarry Smith DM *subdomain_dm = NULL; 3408f3b4b4dSDmitry Karpeev char **subdomain_names = NULL; 341f771a274SFande Kong PetscInt *numbering; 3428f3b4b4dSDmitry Karpeev 343b1a0a8a3SJed Brown 344b1a0a8a3SJed Brown PetscFunctionBegin; 345ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr); 346ce94432eSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);CHKERRQ(ierr); 347b1a0a8a3SJed Brown if (!pc->setupcalled) { 348ea91fabdSFande Kong /* use a hierarchical partitioning */ 349ea91fabdSFande Kong if(osm->hierarchicalpartitioning){ 350ea91fabdSFande Kong ierr = PCGASMSetHierarchicalPartitioning(pc);CHKERRQ(ierr); 351ea91fabdSFande Kong } 352b1a0a8a3SJed Brown if (!osm->type_set) { 353b1a0a8a3SJed Brown ierr = MatIsSymmetricKnown(pc->pmat,&symset,&flg);CHKERRQ(ierr); 3542fa5cd67SKarl Rupp if (symset && flg) osm->type = PC_GASM_BASIC; 355b1a0a8a3SJed Brown } 356b1a0a8a3SJed Brown 3578f3b4b4dSDmitry Karpeev if (osm->n == PETSC_DETERMINE) { 3588f3b4b4dSDmitry Karpeev if (osm->N != PETSC_DETERMINE) { 3598f3b4b4dSDmitry Karpeev /* No local subdomains given, but the desired number of total subdomains is known, so construct them accordingly. */ 3608f3b4b4dSDmitry Karpeev ierr = PCGASMCreateSubdomains(pc->pmat,osm->N,&osm->n,&osm->iis);CHKERRQ(ierr); 3618f3b4b4dSDmitry Karpeev } else if (osm->dm_subdomains && pc->dm) { 3628f3b4b4dSDmitry Karpeev /* try pc->dm next, if allowed */ 3638f3b4b4dSDmitry Karpeev PetscInt d; 364a35b7b57SDmitry Karpeev IS *inner_subdomain_is, *outer_subdomain_is; 365a35b7b57SDmitry Karpeev ierr = DMCreateDomainDecomposition(pc->dm, &num_subdomains, &subdomain_names, &inner_subdomain_is, &outer_subdomain_is, &subdomain_dm);CHKERRQ(ierr); 366a35b7b57SDmitry Karpeev if (num_subdomains) { 367a35b7b57SDmitry Karpeev ierr = PCGASMSetSubdomains(pc, num_subdomains, inner_subdomain_is, outer_subdomain_is);CHKERRQ(ierr); 36869ca1f37SDmitry Karpeev } 369a35b7b57SDmitry Karpeev for (d = 0; d < num_subdomains; ++d) { 370a35b7b57SDmitry Karpeev if (inner_subdomain_is) {ierr = ISDestroy(&inner_subdomain_is[d]);CHKERRQ(ierr);} 371a35b7b57SDmitry Karpeev if (outer_subdomain_is) {ierr = ISDestroy(&outer_subdomain_is[d]);CHKERRQ(ierr);} 372fdc48646SDmitry Karpeev } 373a35b7b57SDmitry Karpeev ierr = PetscFree(inner_subdomain_is);CHKERRQ(ierr); 374a35b7b57SDmitry Karpeev ierr = PetscFree(outer_subdomain_is);CHKERRQ(ierr); 3758f3b4b4dSDmitry Karpeev } else { 3768f3b4b4dSDmitry Karpeev /* still no subdomains; use one per processor */ 37744fe18e5SDmitry Karpeev osm->nmax = osm->n = 1; 378ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr); 379f746d493SDmitry Karpeev osm->N = size; 3808f3b4b4dSDmitry Karpeev ierr = PCGASMCreateLocalSubdomains(pc->pmat,osm->n,&osm->iis);CHKERRQ(ierr); 381fdc48646SDmitry Karpeev } 38206b43e7eSDmitry Karpeev } 383a35b7b57SDmitry Karpeev if (!osm->iis) { 384a35b7b57SDmitry Karpeev /* 3858f3b4b4dSDmitry Karpeev osm->n was set in PCGASMSetSubdomains(), but the actual subdomains have not been supplied. 3868f3b4b4dSDmitry Karpeev We create the requisite number of local inner subdomains and then expand them into 3878f3b4b4dSDmitry Karpeev out subdomains, if necessary. 388a35b7b57SDmitry Karpeev */ 3898f3b4b4dSDmitry Karpeev ierr = PCGASMCreateLocalSubdomains(pc->pmat,osm->n,&osm->iis);CHKERRQ(ierr); 390a35b7b57SDmitry Karpeev } 3918f3b4b4dSDmitry Karpeev if (!osm->ois) { 3928f3b4b4dSDmitry Karpeev /* 3938f3b4b4dSDmitry Karpeev Initially make outer subdomains the same as inner subdomains. If nonzero additional overlap 3948f3b4b4dSDmitry Karpeev has been requested, copy the inner subdomains over so they can be modified. 3958f3b4b4dSDmitry Karpeev */ 3968f3b4b4dSDmitry Karpeev ierr = PetscMalloc1(osm->n,&osm->ois);CHKERRQ(ierr); 3978f3b4b4dSDmitry Karpeev for (i=0; i<osm->n; ++i) { 398ea91fabdSFande Kong if (osm->overlap > 0 && osm->N>1) { /* With positive overlap, osm->iis[i] will be modified */ 3998f3b4b4dSDmitry Karpeev ierr = ISDuplicate(osm->iis[i],(osm->ois)+i);CHKERRQ(ierr); 4008f3b4b4dSDmitry Karpeev ierr = ISCopy(osm->iis[i],osm->ois[i]);CHKERRQ(ierr); 4018f3b4b4dSDmitry Karpeev } else { 4028f3b4b4dSDmitry Karpeev ierr = PetscObjectReference((PetscObject)((osm->iis)[i]));CHKERRQ(ierr); 4038f3b4b4dSDmitry Karpeev osm->ois[i] = osm->iis[i]; 4048f3b4b4dSDmitry Karpeev } 4058f3b4b4dSDmitry Karpeev } 406ea91fabdSFande Kong if (osm->overlap>0 && osm->N>1) { 4078f3b4b4dSDmitry Karpeev /* Extend the "overlapping" regions by a number of steps */ 408d21f2a47SFande Kong ierr = MatIncreaseOverlapSplit(pc->pmat,osm->n,osm->ois,osm->overlap);CHKERRQ(ierr); 4098f3b4b4dSDmitry Karpeev } 410b1a0a8a3SJed Brown } 4116a4f0f73SDmitry Karpeev 4128f3b4b4dSDmitry Karpeev /* Now the subdomains are defined. Determine their global and max local numbers, if necessary. */ 4138f3b4b4dSDmitry Karpeev if (osm->nmax == PETSC_DETERMINE) { 4148f3b4b4dSDmitry Karpeev PetscMPIInt inwork,outwork; 4158f3b4b4dSDmitry Karpeev /* determine global number of subdomains and the max number of local subdomains */ 4168f3b4b4dSDmitry Karpeev inwork = osm->n; 417b2566f29SBarry Smith ierr = MPIU_Allreduce(&inwork,&outwork,1,MPI_INT,MPI_MAX,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr); 4188f3b4b4dSDmitry Karpeev osm->nmax = outwork; 4198f3b4b4dSDmitry Karpeev } 4208f3b4b4dSDmitry Karpeev if (osm->N == PETSC_DETERMINE) { 4218f3b4b4dSDmitry Karpeev /* Determine the number of globally-distinct subdomains and compute a global numbering for them. */ 4228f3b4b4dSDmitry Karpeev ierr = PetscObjectsListGetGlobalNumbering(PetscObjectComm((PetscObject)pc),osm->n,(PetscObject*)osm->ois,&osm->N,NULL);CHKERRQ(ierr); 4238f3b4b4dSDmitry Karpeev } 4248f3b4b4dSDmitry Karpeev 425b1a0a8a3SJed Brown 426b1a0a8a3SJed Brown if (osm->sort_indices) { 427f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 4286a4f0f73SDmitry Karpeev ierr = ISSort(osm->ois[i]);CHKERRQ(ierr); 4296a4f0f73SDmitry Karpeev ierr = ISSort(osm->iis[i]);CHKERRQ(ierr); 430b1a0a8a3SJed Brown } 431b1a0a8a3SJed Brown } 4328f3b4b4dSDmitry Karpeev ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 4338f3b4b4dSDmitry Karpeev ierr = PCGASMPrintSubdomains(pc);CHKERRQ(ierr); 4348f3b4b4dSDmitry Karpeev 4356a4f0f73SDmitry Karpeev /* 4366a4f0f73SDmitry Karpeev Merge the ISs, create merged vectors and restrictions. 4376a4f0f73SDmitry Karpeev */ 4386a4f0f73SDmitry Karpeev /* Merge outer subdomain ISs and construct a restriction onto the disjoint union of local outer subdomains. */ 4396a4f0f73SDmitry Karpeev on = 0; 440f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 4416a4f0f73SDmitry Karpeev ierr = ISGetLocalSize(osm->ois[i],&oni);CHKERRQ(ierr); 4426a4f0f73SDmitry Karpeev on += oni; 443f746d493SDmitry Karpeev } 444785e854fSJed Brown ierr = PetscMalloc1(on, &oidx);CHKERRQ(ierr); 4456a4f0f73SDmitry Karpeev on = 0; 446930d09c1SFande Kong /* Merge local indices together */ 447f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 4486a4f0f73SDmitry Karpeev ierr = ISGetLocalSize(osm->ois[i],&oni);CHKERRQ(ierr); 4496a4f0f73SDmitry Karpeev ierr = ISGetIndices(osm->ois[i],&oidxi);CHKERRQ(ierr); 4506a4f0f73SDmitry Karpeev ierr = PetscMemcpy(oidx+on,oidxi,sizeof(PetscInt)*oni);CHKERRQ(ierr); 4516a4f0f73SDmitry Karpeev ierr = ISRestoreIndices(osm->ois[i],&oidxi);CHKERRQ(ierr); 4526a4f0f73SDmitry Karpeev on += oni; 453f746d493SDmitry Karpeev } 4546a4f0f73SDmitry Karpeev ierr = ISCreateGeneral(((PetscObject)(pc))->comm,on,oidx,PETSC_OWN_POINTER,&gois);CHKERRQ(ierr); 455ea91fabdSFande Kong nTotalInnerIndices = 0; 456ea91fabdSFande Kong for(i=0; i<osm->n; i++){ 457ea91fabdSFande Kong ierr = ISGetLocalSize(osm->iis[i],&nInnerIndices);CHKERRQ(ierr); 458ea91fabdSFande Kong nTotalInnerIndices += nInnerIndices; 459ea91fabdSFande Kong } 460ea91fabdSFande Kong ierr = VecCreateMPI(((PetscObject)(pc))->comm,nTotalInnerIndices,PETSC_DETERMINE,&x);CHKERRQ(ierr); 461ea91fabdSFande Kong ierr = VecDuplicate(x,&y);CHKERRQ(ierr); 462ea91fabdSFande Kong 463ce94432eSBarry Smith ierr = VecCreateMPI(PetscObjectComm((PetscObject)pc),on,PETSC_DECIDE,&osm->gx);CHKERRQ(ierr); 464f746d493SDmitry Karpeev ierr = VecDuplicate(osm->gx,&osm->gy);CHKERRQ(ierr); 465930d09c1SFande Kong ierr = VecGetOwnershipRange(osm->gx, &gostart, NULL);CHKERRQ(ierr); 466930d09c1SFande Kong ierr = ISCreateStride(PetscObjectComm((PetscObject)pc),on,gostart,1, &goid);CHKERRQ(ierr); 467930d09c1SFande Kong /* gois might indices not on local */ 4689448b7f1SJunchao Zhang ierr = VecScatterCreate(x,gois,osm->gx,goid, &(osm->gorestriction));CHKERRQ(ierr); 469f771a274SFande Kong ierr = PetscMalloc1(osm->n,&numbering);CHKERRQ(ierr); 470f771a274SFande Kong ierr = PetscObjectsListGetGlobalNumbering(PetscObjectComm((PetscObject)pc),osm->n,(PetscObject*)osm->ois,NULL,numbering);CHKERRQ(ierr); 4716a4f0f73SDmitry Karpeev ierr = VecDestroy(&x);CHKERRQ(ierr); 4727105d80bSDmitry Karpeev ierr = ISDestroy(&gois);CHKERRQ(ierr); 4732fa5cd67SKarl Rupp 4746a4f0f73SDmitry Karpeev /* Merge inner subdomain ISs and construct a restriction onto the disjoint union of local inner subdomains. */ 4752fa5cd67SKarl Rupp { 4762fa5cd67SKarl Rupp PetscInt ini; /* Number of indices the i-th a local inner subdomain. */ 4776a4f0f73SDmitry Karpeev PetscInt in; /* Number of indices in the disjoint uniont of local inner subdomains. */ 4786a4f0f73SDmitry Karpeev PetscInt *iidx; /* Global indices in the merged local inner subdomain. */ 4796a4f0f73SDmitry Karpeev PetscInt *ioidx; /* Global indices of the disjoint union of inner subdomains within the disjoint union of outer subdomains. */ 4806a4f0f73SDmitry Karpeev IS giis; /* IS for the disjoint union of inner subdomains. */ 4816a4f0f73SDmitry Karpeev IS giois; /* IS for the disjoint union of inner subdomains within the disjoint union of outer subdomains. */ 482f771a274SFande Kong PetscScalar *array; 483f771a274SFande Kong const PetscInt *indices; 484f771a274SFande Kong PetscInt k; 4856a4f0f73SDmitry Karpeev on = 0; 486f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 4876a4f0f73SDmitry Karpeev ierr = ISGetLocalSize(osm->ois[i],&oni);CHKERRQ(ierr); 4886a4f0f73SDmitry Karpeev on += oni; 489f746d493SDmitry Karpeev } 490f771a274SFande Kong ierr = PetscMalloc1(on, &iidx);CHKERRQ(ierr); 491f771a274SFande Kong ierr = PetscMalloc1(on, &ioidx);CHKERRQ(ierr); 492f771a274SFande Kong ierr = VecGetArray(y,&array);CHKERRQ(ierr); 493e12b4729SFande Kong /* set communicator id to determine where overlap is */ 494f771a274SFande Kong in = 0; 495f771a274SFande Kong for (i=0; i<osm->n; i++) { 496f771a274SFande Kong ierr = ISGetLocalSize(osm->iis[i],&ini);CHKERRQ(ierr); 497f771a274SFande Kong for (k = 0; k < ini; ++k){ 498f771a274SFande Kong array[in+k] = numbering[i]; 499f771a274SFande Kong } 500f771a274SFande Kong in += ini; 501f771a274SFande Kong } 502f771a274SFande Kong ierr = VecRestoreArray(y,&array);CHKERRQ(ierr); 503f771a274SFande Kong ierr = VecScatterBegin(osm->gorestriction,y,osm->gy,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 504f771a274SFande Kong ierr = VecScatterEnd(osm->gorestriction,y,osm->gy,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 505f771a274SFande Kong ierr = VecGetOwnershipRange(osm->gy,&gostart, NULL);CHKERRQ(ierr); 506f771a274SFande Kong ierr = VecGetArray(osm->gy,&array);CHKERRQ(ierr); 507f771a274SFande Kong on = 0; 508f771a274SFande Kong in = 0; 509f771a274SFande Kong for (i=0; i<osm->n; i++) { 510f771a274SFande Kong ierr = ISGetLocalSize(osm->ois[i],&oni);CHKERRQ(ierr); 511f771a274SFande Kong ierr = ISGetIndices(osm->ois[i],&indices);CHKERRQ(ierr); 512f771a274SFande Kong for (k=0; k<oni; k++) { 513e12b4729SFande Kong /* skip overlapping indices to get inner domain */ 51443081b6cSDmitry Karpeev if(PetscRealPart(array[on+k]) != numbering[i]) continue; 515f771a274SFande Kong iidx[in] = indices[k]; 516f771a274SFande Kong ioidx[in++] = gostart+on+k; 517f771a274SFande Kong } 518f771a274SFande Kong ierr = ISRestoreIndices(osm->ois[i], &indices);CHKERRQ(ierr); 519f771a274SFande Kong on += oni; 520f771a274SFande Kong } 521f771a274SFande Kong ierr = VecRestoreArray(osm->gy,&array);CHKERRQ(ierr); 522ce94432eSBarry Smith ierr = ISCreateGeneral(PetscObjectComm((PetscObject)pc),in,iidx,PETSC_OWN_POINTER,&giis);CHKERRQ(ierr); 523ce94432eSBarry Smith ierr = ISCreateGeneral(PetscObjectComm((PetscObject)pc),in,ioidx,PETSC_OWN_POINTER,&giois);CHKERRQ(ierr); 5249448b7f1SJunchao Zhang ierr = VecScatterCreate(y,giis,osm->gy,giois,&osm->girestriction);CHKERRQ(ierr); 5256a4f0f73SDmitry Karpeev ierr = VecDestroy(&y);CHKERRQ(ierr); 5266a4f0f73SDmitry Karpeev ierr = ISDestroy(&giis);CHKERRQ(ierr); 5276a4f0f73SDmitry Karpeev ierr = ISDestroy(&giois);CHKERRQ(ierr); 528b1a0a8a3SJed Brown } 5296a4f0f73SDmitry Karpeev ierr = ISDestroy(&goid);CHKERRQ(ierr); 530f771a274SFande Kong ierr = PetscFree(numbering);CHKERRQ(ierr); 5312fa5cd67SKarl Rupp 532f746d493SDmitry Karpeev /* Create the subdomain work vectors. */ 533785e854fSJed Brown ierr = PetscMalloc1(osm->n,&osm->x);CHKERRQ(ierr); 534785e854fSJed Brown ierr = PetscMalloc1(osm->n,&osm->y);CHKERRQ(ierr); 535f746d493SDmitry Karpeev ierr = VecGetArray(osm->gx, &gxarray);CHKERRQ(ierr); 536f746d493SDmitry Karpeev ierr = VecGetArray(osm->gy, &gyarray);CHKERRQ(ierr); 5376a4f0f73SDmitry Karpeev for (i=0, on=0; i<osm->n; ++i, on += oni) { 5386a4f0f73SDmitry Karpeev PetscInt oNi; 5396a4f0f73SDmitry Karpeev ierr = ISGetLocalSize(osm->ois[i],&oni);CHKERRQ(ierr); 540930d09c1SFande Kong /* on a sub communicator */ 5416a4f0f73SDmitry Karpeev ierr = ISGetSize(osm->ois[i],&oNi);CHKERRQ(ierr); 5426a4f0f73SDmitry Karpeev ierr = VecCreateMPIWithArray(((PetscObject)(osm->ois[i]))->comm,1,oni,oNi,gxarray+on,&osm->x[i]);CHKERRQ(ierr); 5436a4f0f73SDmitry Karpeev ierr = VecCreateMPIWithArray(((PetscObject)(osm->ois[i]))->comm,1,oni,oNi,gyarray+on,&osm->y[i]);CHKERRQ(ierr); 544b1a0a8a3SJed Brown } 545f746d493SDmitry Karpeev ierr = VecRestoreArray(osm->gx, &gxarray);CHKERRQ(ierr); 546f746d493SDmitry Karpeev ierr = VecRestoreArray(osm->gy, &gyarray);CHKERRQ(ierr); 5478f3b4b4dSDmitry Karpeev /* Create the subdomain solvers */ 548785e854fSJed Brown ierr = PetscMalloc1(osm->n,&osm->ksp);CHKERRQ(ierr); 5498f3b4b4dSDmitry Karpeev for (i=0; i<osm->n; i++) { 5508f3b4b4dSDmitry Karpeev char subprefix[PETSC_MAX_PATH_LEN+1]; 5516a4f0f73SDmitry Karpeev ierr = KSPCreate(((PetscObject)(osm->ois[i]))->comm,&ksp);CHKERRQ(ierr); 55266b14c46SBarry Smith ierr = KSPSetErrorIfNotConverged(ksp,pc->erroriffailure);CHKERRQ(ierr); 5533bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ksp);CHKERRQ(ierr); 554b1a0a8a3SJed Brown ierr = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,1);CHKERRQ(ierr); 555b1a0a8a3SJed Brown ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr); 5568f3b4b4dSDmitry Karpeev ierr = KSPGetPC(ksp,&subpc);CHKERRQ(ierr); /* Why do we need this here? */ 5578f3b4b4dSDmitry Karpeev if (subdomain_dm) { 5588f3b4b4dSDmitry Karpeev ierr = KSPSetDM(ksp,subdomain_dm[i]);CHKERRQ(ierr); 5598f3b4b4dSDmitry Karpeev ierr = DMDestroy(subdomain_dm+i);CHKERRQ(ierr); 5608f3b4b4dSDmitry Karpeev } 561b1a0a8a3SJed Brown ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 562b1a0a8a3SJed Brown ierr = KSPSetOptionsPrefix(ksp,prefix);CHKERRQ(ierr); 5638f3b4b4dSDmitry Karpeev if (subdomain_names && subdomain_names[i]) { 5648f3b4b4dSDmitry Karpeev ierr = PetscSNPrintf(subprefix,PETSC_MAX_PATH_LEN,"sub_%s_",subdomain_names[i]);CHKERRQ(ierr); 5658f3b4b4dSDmitry Karpeev ierr = KSPAppendOptionsPrefix(ksp,subprefix);CHKERRQ(ierr); 5668f3b4b4dSDmitry Karpeev ierr = PetscFree(subdomain_names[i]);CHKERRQ(ierr); 5678f3b4b4dSDmitry Karpeev } 568b1a0a8a3SJed Brown ierr = KSPAppendOptionsPrefix(ksp,"sub_");CHKERRQ(ierr); 569b1a0a8a3SJed Brown osm->ksp[i] = ksp; 570b1a0a8a3SJed Brown } 5718f3b4b4dSDmitry Karpeev ierr = PetscFree(subdomain_dm);CHKERRQ(ierr); 5728f3b4b4dSDmitry Karpeev ierr = PetscFree(subdomain_names);CHKERRQ(ierr); 573b1a0a8a3SJed Brown scall = MAT_INITIAL_MATRIX; 574b1a0a8a3SJed Brown 5758f3b4b4dSDmitry Karpeev } else { /* if (pc->setupcalled) */ 576b1a0a8a3SJed Brown /* 5778f3b4b4dSDmitry Karpeev Destroy the submatrices from the previous iteration 578b1a0a8a3SJed Brown */ 579b1a0a8a3SJed Brown if (pc->flag == DIFFERENT_NONZERO_PATTERN) { 580f746d493SDmitry Karpeev ierr = MatDestroyMatrices(osm->n,&osm->pmat);CHKERRQ(ierr); 581b1a0a8a3SJed Brown scall = MAT_INITIAL_MATRIX; 582b1a0a8a3SJed Brown } 583ea91fabdSFande Kong if(osm->permutationIS){ 5847dae84e0SHong Zhang ierr = MatCreateSubMatrix(pc->pmat,osm->permutationIS,osm->permutationIS,scall,&osm->permutationP);CHKERRQ(ierr); 585ea91fabdSFande Kong ierr = PetscObjectReference((PetscObject)osm->permutationP);CHKERRQ(ierr); 586ea91fabdSFande Kong osm->pcmat = pc->pmat; 587ea91fabdSFande Kong pc->pmat = osm->permutationP; 588b1a0a8a3SJed Brown } 589b1a0a8a3SJed Brown 590ea91fabdSFande Kong } 591ea91fabdSFande Kong 592ea91fabdSFande Kong 593b1a0a8a3SJed Brown /* 594f746d493SDmitry Karpeev Extract out the submatrices. 595b1a0a8a3SJed Brown */ 59681a5c4aaSDmitry Karpeev if (size > 1) { 5977dae84e0SHong Zhang ierr = MatCreateSubMatricesMPI(pc->pmat,osm->n,osm->ois,osm->ois,scall,&osm->pmat);CHKERRQ(ierr); 5982fa5cd67SKarl Rupp } else { 5997dae84e0SHong Zhang ierr = MatCreateSubMatrices(pc->pmat,osm->n,osm->ois,osm->ois,scall,&osm->pmat);CHKERRQ(ierr); 60081a5c4aaSDmitry Karpeev } 601b1a0a8a3SJed Brown if (scall == MAT_INITIAL_MATRIX) { 602b1a0a8a3SJed Brown ierr = PetscObjectGetOptionsPrefix((PetscObject)pc->pmat,&pprefix);CHKERRQ(ierr); 603f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 6043bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)osm->pmat[i]);CHKERRQ(ierr); 605b1a0a8a3SJed Brown ierr = PetscObjectSetOptionsPrefix((PetscObject)osm->pmat[i],pprefix);CHKERRQ(ierr); 606b1a0a8a3SJed Brown } 607b1a0a8a3SJed Brown } 608b1a0a8a3SJed Brown 609b1a0a8a3SJed Brown /* Return control to the user so that the submatrices can be modified (e.g., to apply 610b1a0a8a3SJed Brown different boundary conditions for the submatrices than for the global problem) */ 6116a4f0f73SDmitry Karpeev ierr = PCModifySubMatrices(pc,osm->n,osm->ois,osm->ois,osm->pmat,pc->modifysubmatricesP);CHKERRQ(ierr); 612b1a0a8a3SJed Brown 613b1a0a8a3SJed Brown /* 6146a4f0f73SDmitry Karpeev Loop over submatrices putting them into local ksps 615b1a0a8a3SJed Brown */ 616f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 61723ee1639SBarry Smith ierr = KSPSetOperators(osm->ksp[i],osm->pmat[i],osm->pmat[i]);CHKERRQ(ierr); 618b1a0a8a3SJed Brown if (!pc->setupcalled) { 619b1a0a8a3SJed Brown ierr = KSPSetFromOptions(osm->ksp[i]);CHKERRQ(ierr); 620b1a0a8a3SJed Brown } 621b1a0a8a3SJed Brown } 622ea91fabdSFande Kong if(osm->pcmat){ 623ea91fabdSFande Kong ierr = MatDestroy(&pc->pmat);CHKERRQ(ierr); 624ea91fabdSFande Kong pc->pmat = osm->pcmat; 625ea91fabdSFande Kong osm->pcmat = 0; 626ea91fabdSFande Kong } 627b1a0a8a3SJed Brown PetscFunctionReturn(0); 628b1a0a8a3SJed Brown } 629b1a0a8a3SJed Brown 630f746d493SDmitry Karpeev static PetscErrorCode PCSetUpOnBlocks_GASM(PC pc) 631b1a0a8a3SJed Brown { 632f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 633b1a0a8a3SJed Brown PetscErrorCode ierr; 634b1a0a8a3SJed Brown PetscInt i; 635b1a0a8a3SJed Brown 636b1a0a8a3SJed Brown PetscFunctionBegin; 637f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 638b1a0a8a3SJed Brown ierr = KSPSetUp(osm->ksp[i]);CHKERRQ(ierr); 639b1a0a8a3SJed Brown } 640b1a0a8a3SJed Brown PetscFunctionReturn(0); 641b1a0a8a3SJed Brown } 642b1a0a8a3SJed Brown 643ea91fabdSFande Kong static PetscErrorCode PCApply_GASM(PC pc,Vec xin,Vec yout) 644b1a0a8a3SJed Brown { 645f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 646b1a0a8a3SJed Brown PetscErrorCode ierr; 647f746d493SDmitry Karpeev PetscInt i; 648ea91fabdSFande Kong Vec x,y; 649b1a0a8a3SJed Brown ScatterMode forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE; 650b1a0a8a3SJed Brown 651b1a0a8a3SJed Brown PetscFunctionBegin; 652ea91fabdSFande Kong if(osm->pctoouter){ 653ea91fabdSFande Kong ierr = VecScatterBegin(osm->pctoouter,xin,osm->pcx,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 654ea91fabdSFande Kong ierr = VecScatterEnd(osm->pctoouter,xin,osm->pcx,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 655ea91fabdSFande Kong x = osm->pcx; 656ea91fabdSFande Kong y = osm->pcy; 657ea91fabdSFande Kong }else{ 658ea91fabdSFande Kong x = xin; 659ea91fabdSFande Kong y = yout; 660ea91fabdSFande Kong } 661b1a0a8a3SJed Brown /* 6626a4f0f73SDmitry Karpeev Support for limiting the restriction or interpolation only to the inner 663b1a0a8a3SJed Brown subdomain values (leaving the other values 0). 664b1a0a8a3SJed Brown */ 665f746d493SDmitry Karpeev if (!(osm->type & PC_GASM_RESTRICT)) { 666b1a0a8a3SJed Brown /* have to zero the work RHS since scatter may leave some slots empty */ 667f746d493SDmitry Karpeev ierr = VecZeroEntries(osm->gx);CHKERRQ(ierr); 6686a4f0f73SDmitry Karpeev ierr = VecScatterBegin(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 6692fa5cd67SKarl Rupp } else { 6706a4f0f73SDmitry Karpeev ierr = VecScatterBegin(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 671b1a0a8a3SJed Brown } 6726a4f0f73SDmitry Karpeev ierr = VecZeroEntries(osm->gy);CHKERRQ(ierr); 6736a4f0f73SDmitry Karpeev if (!(osm->type & PC_GASM_RESTRICT)) { 6746a4f0f73SDmitry Karpeev ierr = VecScatterEnd(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 6752fa5cd67SKarl Rupp } else { 6766a4f0f73SDmitry Karpeev ierr = VecScatterEnd(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 6776a4f0f73SDmitry Karpeev } 67886b96d36SDmitry Karpeev /* do the subdomain solves */ 67986b96d36SDmitry Karpeev for (i=0; i<osm->n; ++i) { 680b1a0a8a3SJed Brown ierr = KSPSolve(osm->ksp[i],osm->x[i],osm->y[i]);CHKERRQ(ierr); 681c0decd05SBarry Smith ierr = KSPCheckSolve(osm->ksp[i],pc,osm->y[i]);CHKERRQ(ierr); 682b1a0a8a3SJed Brown } 683930d09c1SFande Kong /* Do we need to zero y ?? */ 6846a4f0f73SDmitry Karpeev ierr = VecZeroEntries(y);CHKERRQ(ierr); 6856a4f0f73SDmitry Karpeev if (!(osm->type & PC_GASM_INTERPOLATE)) { 6866a4f0f73SDmitry Karpeev ierr = VecScatterBegin(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); 687ea91fabdSFande Kong ierr = VecScatterEnd(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); 6882fa5cd67SKarl Rupp } else { 6896a4f0f73SDmitry Karpeev ierr = VecScatterBegin(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); 690ea91fabdSFande Kong ierr = VecScatterEnd(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); 6916a4f0f73SDmitry Karpeev } 692ea91fabdSFande Kong if(osm->pctoouter){ 693ea91fabdSFande Kong ierr = VecScatterBegin(osm->pctoouter,y,yout,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 694ea91fabdSFande Kong ierr = VecScatterEnd(osm->pctoouter,y,yout,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 695ea91fabdSFande Kong } 696ea91fabdSFande Kong PetscFunctionReturn(0); 697b1a0a8a3SJed Brown } 698b1a0a8a3SJed Brown 699ea91fabdSFande Kong static PetscErrorCode PCApplyTranspose_GASM(PC pc,Vec xin,Vec yout) 700b1a0a8a3SJed Brown { 701f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 702b1a0a8a3SJed Brown PetscErrorCode ierr; 703f746d493SDmitry Karpeev PetscInt i; 704ea91fabdSFande Kong Vec x,y; 705b1a0a8a3SJed Brown ScatterMode forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE; 706b1a0a8a3SJed Brown 707b1a0a8a3SJed Brown PetscFunctionBegin; 708ea91fabdSFande Kong if(osm->pctoouter){ 709ea91fabdSFande Kong ierr = VecScatterBegin(osm->pctoouter,xin,osm->pcx,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 710ea91fabdSFande Kong ierr = VecScatterEnd(osm->pctoouter,xin,osm->pcx,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 711ea91fabdSFande Kong x = osm->pcx; 712ea91fabdSFande Kong y = osm->pcy; 713ea91fabdSFande Kong }else{ 714ea91fabdSFande Kong x = xin; 715ea91fabdSFande Kong y = yout; 716ea91fabdSFande Kong } 717b1a0a8a3SJed Brown /* 718b1a0a8a3SJed Brown Support for limiting the restriction or interpolation to only local 719b1a0a8a3SJed Brown subdomain values (leaving the other values 0). 720b1a0a8a3SJed Brown 721f746d493SDmitry Karpeev Note: these are reversed from the PCApply_GASM() because we are applying the 722b1a0a8a3SJed Brown transpose of the three terms 723b1a0a8a3SJed Brown */ 724f746d493SDmitry Karpeev if (!(osm->type & PC_GASM_INTERPOLATE)) { 725b1a0a8a3SJed Brown /* have to zero the work RHS since scatter may leave some slots empty */ 726f746d493SDmitry Karpeev ierr = VecZeroEntries(osm->gx);CHKERRQ(ierr); 7276a4f0f73SDmitry Karpeev ierr = VecScatterBegin(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 7282fa5cd67SKarl Rupp } else { 7296a4f0f73SDmitry Karpeev ierr = VecScatterBegin(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 730b1a0a8a3SJed Brown } 7316a4f0f73SDmitry Karpeev ierr = VecZeroEntries(osm->gy);CHKERRQ(ierr); 7326a4f0f73SDmitry Karpeev if (!(osm->type & PC_GASM_INTERPOLATE)) { 7336a4f0f73SDmitry Karpeev ierr = VecScatterEnd(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 7342fa5cd67SKarl Rupp } else { 7356a4f0f73SDmitry Karpeev ierr = VecScatterEnd(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 7366a4f0f73SDmitry Karpeev } 737b1a0a8a3SJed Brown /* do the local solves */ 738ab3e923fSDmitry 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. */ 739b1a0a8a3SJed Brown ierr = KSPSolveTranspose(osm->ksp[i],osm->x[i],osm->y[i]);CHKERRQ(ierr); 740c0decd05SBarry Smith ierr = KSPCheckSolve(osm->ksp[i],pc,osm->y[i]);CHKERRQ(ierr); 741b1a0a8a3SJed Brown } 7426a4f0f73SDmitry Karpeev ierr = VecZeroEntries(y);CHKERRQ(ierr); 7436a4f0f73SDmitry Karpeev if (!(osm->type & PC_GASM_RESTRICT)) { 7446a4f0f73SDmitry Karpeev ierr = VecScatterBegin(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); 7456a4f0f73SDmitry Karpeev ierr = VecScatterEnd(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); 7462fa5cd67SKarl Rupp } else { 7476a4f0f73SDmitry Karpeev ierr = VecScatterBegin(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); 7486a4f0f73SDmitry Karpeev ierr = VecScatterEnd(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); 7496a4f0f73SDmitry Karpeev } 750ea91fabdSFande Kong if(osm->pctoouter){ 751ea91fabdSFande Kong ierr = VecScatterBegin(osm->pctoouter,y,yout,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 752ea91fabdSFande Kong ierr = VecScatterEnd(osm->pctoouter,y,yout,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 753ea91fabdSFande Kong } 754b1a0a8a3SJed Brown PetscFunctionReturn(0); 755b1a0a8a3SJed Brown } 756b1a0a8a3SJed Brown 757a06653b4SBarry Smith static PetscErrorCode PCReset_GASM(PC pc) 758b1a0a8a3SJed Brown { 759f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 760b1a0a8a3SJed Brown PetscErrorCode ierr; 761b1a0a8a3SJed Brown PetscInt i; 762b1a0a8a3SJed Brown 763b1a0a8a3SJed Brown PetscFunctionBegin; 764b1a0a8a3SJed Brown if (osm->ksp) { 765f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 766a06653b4SBarry Smith ierr = KSPReset(osm->ksp[i]);CHKERRQ(ierr); 767b1a0a8a3SJed Brown } 768b1a0a8a3SJed Brown } 769b1a0a8a3SJed Brown if (osm->pmat) { 770f746d493SDmitry Karpeev if (osm->n > 0) { 771df750dc8SHong Zhang PetscMPIInt size; 772df750dc8SHong Zhang ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr); 773df750dc8SHong Zhang if (size > 1) { 7747dae84e0SHong Zhang /* osm->pmat is created by MatCreateSubMatricesMPI(), cannot use MatDestroySubMatrices() */ 775ab3e923fSDmitry Karpeev ierr = MatDestroyMatrices(osm->n,&osm->pmat);CHKERRQ(ierr); 776df750dc8SHong Zhang } else { 777df750dc8SHong Zhang ierr = MatDestroySubMatrices(osm->n,&osm->pmat);CHKERRQ(ierr); 778df750dc8SHong Zhang } 779b1a0a8a3SJed Brown } 780b1a0a8a3SJed Brown } 781d34fcf5fSBarry Smith if (osm->x) { 782f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 783fcfd50ebSBarry Smith ierr = VecDestroy(&osm->x[i]);CHKERRQ(ierr); 784fcfd50ebSBarry Smith ierr = VecDestroy(&osm->y[i]);CHKERRQ(ierr); 785b1a0a8a3SJed Brown } 786d34fcf5fSBarry Smith } 7876bf464f9SBarry Smith ierr = VecDestroy(&osm->gx);CHKERRQ(ierr); 7886bf464f9SBarry Smith ierr = VecDestroy(&osm->gy);CHKERRQ(ierr); 789ab3e923fSDmitry Karpeev 7906a4f0f73SDmitry Karpeev ierr = VecScatterDestroy(&osm->gorestriction);CHKERRQ(ierr); 7916a4f0f73SDmitry Karpeev ierr = VecScatterDestroy(&osm->girestriction);CHKERRQ(ierr); 7928f3b4b4dSDmitry Karpeev if (!osm->user_subdomains) { 7932c112581SDmitry Karpeev ierr = PCGASMDestroySubdomains(osm->n,&osm->ois,&osm->iis);CHKERRQ(ierr); 7948f3b4b4dSDmitry Karpeev osm->N = PETSC_DETERMINE; 7958f3b4b4dSDmitry Karpeev osm->nmax = PETSC_DETERMINE; 7968f3b4b4dSDmitry Karpeev } 797ea91fabdSFande Kong if(osm->pctoouter){ 798ea91fabdSFande Kong ierr = VecScatterDestroy(&(osm->pctoouter));CHKERRQ(ierr); 799ea91fabdSFande Kong } 800ea91fabdSFande Kong if(osm->permutationIS){ 801ea91fabdSFande Kong ierr = ISDestroy(&(osm->permutationIS));CHKERRQ(ierr); 802ea91fabdSFande Kong } 803ea91fabdSFande Kong if(osm->pcx){ 804ea91fabdSFande Kong ierr = VecDestroy(&(osm->pcx));CHKERRQ(ierr); 805ea91fabdSFande Kong } 806ea91fabdSFande Kong if(osm->pcy){ 807ea91fabdSFande Kong ierr = VecDestroy(&(osm->pcy));CHKERRQ(ierr); 808ea91fabdSFande Kong } 809ea91fabdSFande Kong if(osm->permutationP){ 810ea91fabdSFande Kong ierr = MatDestroy(&(osm->permutationP));CHKERRQ(ierr); 811ea91fabdSFande Kong } 812ea91fabdSFande Kong if(osm->pcmat){ 813ea91fabdSFande Kong ierr = MatDestroy(&osm->pcmat);CHKERRQ(ierr); 814ea91fabdSFande Kong } 815a06653b4SBarry Smith PetscFunctionReturn(0); 816a06653b4SBarry Smith } 817a06653b4SBarry Smith 818a06653b4SBarry Smith static PetscErrorCode PCDestroy_GASM(PC pc) 819a06653b4SBarry Smith { 820a06653b4SBarry Smith PC_GASM *osm = (PC_GASM*)pc->data; 821a06653b4SBarry Smith PetscErrorCode ierr; 822a06653b4SBarry Smith PetscInt i; 823a06653b4SBarry Smith 824a06653b4SBarry Smith PetscFunctionBegin; 825a06653b4SBarry Smith ierr = PCReset_GASM(pc);CHKERRQ(ierr); 826135757f6SDmitry Karpeev /* PCReset will not destroy subdomains, if user_subdomains is true. */ 827135757f6SDmitry Karpeev ierr = PCGASMDestroySubdomains(osm->n,&osm->ois,&osm->iis);CHKERRQ(ierr); 828a06653b4SBarry Smith if (osm->ksp) { 829a06653b4SBarry Smith for (i=0; i<osm->n; i++) { 8306bf464f9SBarry Smith ierr = KSPDestroy(&osm->ksp[i]);CHKERRQ(ierr); 831a06653b4SBarry Smith } 832a06653b4SBarry Smith ierr = PetscFree(osm->ksp);CHKERRQ(ierr); 833a06653b4SBarry Smith } 834a06653b4SBarry Smith ierr = PetscFree(osm->x);CHKERRQ(ierr); 835a06653b4SBarry Smith ierr = PetscFree(osm->y);CHKERRQ(ierr); 836c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 837b1a0a8a3SJed Brown PetscFunctionReturn(0); 838b1a0a8a3SJed Brown } 839b1a0a8a3SJed Brown 8404416b707SBarry Smith static PetscErrorCode PCSetFromOptions_GASM(PetscOptionItems *PetscOptionsObject,PC pc) 841a6dfd86eSKarl Rupp { 842f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 843b1a0a8a3SJed Brown PetscErrorCode ierr; 844b1a0a8a3SJed Brown PetscInt blocks,ovl; 845b1a0a8a3SJed Brown PetscBool symset,flg; 846f746d493SDmitry Karpeev PCGASMType gasmtype; 847b1a0a8a3SJed Brown 848b1a0a8a3SJed Brown PetscFunctionBegin; 849b1a0a8a3SJed Brown /* set the type to symmetric if matrix is symmetric */ 850b1a0a8a3SJed Brown if (!osm->type_set && pc->pmat) { 851b1a0a8a3SJed Brown ierr = MatIsSymmetricKnown(pc->pmat,&symset,&flg);CHKERRQ(ierr); 8522fa5cd67SKarl Rupp if (symset && flg) osm->type = PC_GASM_BASIC; 853b1a0a8a3SJed Brown } 854e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"Generalized additive Schwarz options");CHKERRQ(ierr); 8558f3b4b4dSDmitry 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); 8568f3b4b4dSDmitry Karpeev ierr = PetscOptionsInt("-pc_gasm_total_subdomains","Total number of subdomains across communicator","PCGASMSetTotalSubdomains",osm->N,&blocks,&flg);CHKERRQ(ierr); 85765db9045SDmitry Karpeev if (flg) { 8588f3b4b4dSDmitry Karpeev ierr = PCGASMSetTotalSubdomains(pc,blocks);CHKERRQ(ierr); 85965db9045SDmitry Karpeev } 86006b43e7eSDmitry Karpeev ierr = PetscOptionsInt("-pc_gasm_overlap","Number of overlapping degrees of freedom","PCGASMSetOverlap",osm->overlap,&ovl,&flg);CHKERRQ(ierr); 86165db9045SDmitry Karpeev if (flg) { 86265db9045SDmitry Karpeev ierr = PCGASMSetOverlap(pc,ovl);CHKERRQ(ierr); 863d709ea83SDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 86465db9045SDmitry Karpeev } 865b1a0a8a3SJed Brown flg = PETSC_FALSE; 866f746d493SDmitry Karpeev ierr = PetscOptionsEnum("-pc_gasm_type","Type of restriction/extension","PCGASMSetType",PCGASMTypes,(PetscEnum)osm->type,(PetscEnum*)&gasmtype,&flg);CHKERRQ(ierr); 867f746d493SDmitry Karpeev if (flg) {ierr = PCGASMSetType(pc,gasmtype);CHKERRQ(ierr);} 868ea91fabdSFande Kong ierr = PetscOptionsBool("-pc_gasm_use_hierachical_partitioning","use hierarchical partitioning",NULL,osm->hierarchicalpartitioning,&osm->hierarchicalpartitioning,&flg);CHKERRQ(ierr); 869b1a0a8a3SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 870b1a0a8a3SJed Brown PetscFunctionReturn(0); 871b1a0a8a3SJed Brown } 872b1a0a8a3SJed Brown 873b1a0a8a3SJed Brown /*------------------------------------------------------------------------------------*/ 874b1a0a8a3SJed Brown 8758f3b4b4dSDmitry Karpeev /*@ 8768f3b4b4dSDmitry Karpeev PCGASMSetTotalSubdomains - sets the total number of subdomains to use across the 8778f3b4b4dSDmitry Karpeev communicator. 8788f3b4b4dSDmitry Karpeev Logically collective on pc 8798f3b4b4dSDmitry Karpeev 8808f3b4b4dSDmitry Karpeev Input Parameters: 8818f3b4b4dSDmitry Karpeev + pc - the preconditioner 8828f3b4b4dSDmitry Karpeev - N - total number of subdomains 8838f3b4b4dSDmitry Karpeev 8848f3b4b4dSDmitry Karpeev 8858f3b4b4dSDmitry Karpeev Level: beginner 8868f3b4b4dSDmitry Karpeev 8878f3b4b4dSDmitry Karpeev .keywords: PC, ASM, DM, set, subdomains, additive Schwarz 8888f3b4b4dSDmitry Karpeev 8898f3b4b4dSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMSetOverlap() 8908f3b4b4dSDmitry Karpeev PCGASMCreateSubdomains2D() 8918f3b4b4dSDmitry Karpeev @*/ 8928f3b4b4dSDmitry Karpeev PetscErrorCode PCGASMSetTotalSubdomains(PC pc,PetscInt N) 8938f3b4b4dSDmitry Karpeev { 8948f3b4b4dSDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 8958f3b4b4dSDmitry Karpeev PetscMPIInt size,rank; 8968f3b4b4dSDmitry Karpeev PetscErrorCode ierr; 8978f3b4b4dSDmitry Karpeev 8988f3b4b4dSDmitry Karpeev PetscFunctionBegin; 8998f3b4b4dSDmitry 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); 9008f3b4b4dSDmitry Karpeev if (pc->setupcalled) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetTotalSubdomains() should be called before calling PCSetUp()."); 9018f3b4b4dSDmitry Karpeev 9022c112581SDmitry Karpeev ierr = PCGASMDestroySubdomains(osm->n,&osm->iis,&osm->ois);CHKERRQ(ierr); 9038f3b4b4dSDmitry Karpeev osm->ois = osm->iis = NULL; 9048f3b4b4dSDmitry Karpeev 9058f3b4b4dSDmitry Karpeev ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr); 9068f3b4b4dSDmitry Karpeev ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);CHKERRQ(ierr); 9078f3b4b4dSDmitry Karpeev osm->N = N; 9088f3b4b4dSDmitry Karpeev osm->n = PETSC_DETERMINE; 9098f3b4b4dSDmitry Karpeev osm->nmax = PETSC_DETERMINE; 9108f3b4b4dSDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 9118f3b4b4dSDmitry Karpeev PetscFunctionReturn(0); 9128f3b4b4dSDmitry Karpeev } 9138f3b4b4dSDmitry Karpeev 914670417b2SFande Kong 915b541e6a4SDmitry Karpeev static PetscErrorCode PCGASMSetSubdomains_GASM(PC pc,PetscInt n,IS iis[],IS ois[]) 916b1a0a8a3SJed Brown { 917f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 918b1a0a8a3SJed Brown PetscErrorCode ierr; 919b1a0a8a3SJed Brown PetscInt i; 920b1a0a8a3SJed Brown 921b1a0a8a3SJed Brown PetscFunctionBegin; 9228f3b4b4dSDmitry Karpeev if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Each process must have 1 or more subdomains, got n = %D",n); 9238f3b4b4dSDmitry Karpeev if (pc->setupcalled) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetSubdomains() should be called before calling PCSetUp()."); 924b1a0a8a3SJed Brown 9252c112581SDmitry Karpeev ierr = PCGASMDestroySubdomains(osm->n,&osm->iis,&osm->ois);CHKERRQ(ierr); 9268f3b4b4dSDmitry Karpeev osm->iis = osm->ois = NULL; 9278f3b4b4dSDmitry Karpeev osm->n = n; 9288f3b4b4dSDmitry Karpeev osm->N = PETSC_DETERMINE; 9298f3b4b4dSDmitry Karpeev osm->nmax = PETSC_DETERMINE; 930a35b7b57SDmitry Karpeev if (ois) { 931785e854fSJed Brown ierr = PetscMalloc1(n,&osm->ois);CHKERRQ(ierr); 9328f3b4b4dSDmitry Karpeev for (i=0; i<n; i++) { 9338f3b4b4dSDmitry Karpeev ierr = PetscObjectReference((PetscObject)ois[i]);CHKERRQ(ierr); 9348f3b4b4dSDmitry Karpeev osm->ois[i] = ois[i]; 9358f3b4b4dSDmitry Karpeev } 9368f3b4b4dSDmitry Karpeev /* 9378f3b4b4dSDmitry Karpeev Since the user set the outer subdomains, even if nontrivial overlap was requested via PCGASMSetOverlap(), 9388f3b4b4dSDmitry Karpeev it will be ignored. To avoid confusion later on (e.g., when viewing the PC), the overlap size is set to -1. 9398f3b4b4dSDmitry Karpeev */ 940b1a0a8a3SJed Brown osm->overlap = -1; 941670417b2SFande Kong /* inner subdomains must be provided */ 942670417b2SFande Kong if (!iis) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"inner indices have to be provided \n"); 943670417b2SFande Kong }/* end if */ 944a35b7b57SDmitry Karpeev if (iis) { 945785e854fSJed Brown ierr = PetscMalloc1(n,&osm->iis);CHKERRQ(ierr); 9468f3b4b4dSDmitry Karpeev for (i=0; i<n; i++) { 9478f3b4b4dSDmitry Karpeev ierr = PetscObjectReference((PetscObject)iis[i]);CHKERRQ(ierr); 9488f3b4b4dSDmitry Karpeev osm->iis[i] = iis[i]; 9498f3b4b4dSDmitry Karpeev } 950a35b7b57SDmitry Karpeev if (!ois) { 951390e1bf2SBarry Smith osm->ois = NULL; 952670417b2SFande Kong /* if user does not provide outer indices, we will create the corresponding outer indices using osm->overlap =1 in PCSetUp_GASM */ 953670417b2SFande Kong } 954670417b2SFande Kong } 955670417b2SFande Kong #if defined(PETSC_USE_DEBUG) 956670417b2SFande Kong { 957670417b2SFande Kong PetscInt j,rstart,rend,*covered,lsize; 958670417b2SFande Kong const PetscInt *indices; 959670417b2SFande Kong /* check if the inner indices cover and only cover the local portion of the preconditioning matrix */ 960670417b2SFande Kong ierr = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr); 961670417b2SFande Kong ierr = PetscCalloc1(rend-rstart,&covered);CHKERRQ(ierr); 962670417b2SFande Kong /* check if the current processor owns indices from others */ 963a35b7b57SDmitry Karpeev for (i=0; i<n; i++) { 964670417b2SFande Kong ierr = ISGetIndices(osm->iis[i],&indices);CHKERRQ(ierr); 965670417b2SFande Kong ierr = ISGetLocalSize(osm->iis[i],&lsize);CHKERRQ(ierr); 966670417b2SFande Kong for (j=0; j<lsize; j++) { 967670417b2SFande 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]); 968670417b2SFande 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]); 969670417b2SFande Kong else covered[indices[j]-rstart] = 1; 970a35b7b57SDmitry Karpeev } 971670417b2SFande Kong ierr = ISRestoreIndices(osm->iis[i],&indices);CHKERRQ(ierr); 9728f3b4b4dSDmitry Karpeev } 973670417b2SFande Kong /* check if we miss any indices */ 974670417b2SFande Kong for (i=rstart; i<rend; i++) { 975670417b2SFande Kong if (!covered[i-rstart]) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"local entity %d was not covered by inner subdomains",i); 976a35b7b57SDmitry Karpeev } 977670417b2SFande Kong ierr = PetscFree(covered);CHKERRQ(ierr); 978a35b7b57SDmitry Karpeev } 979670417b2SFande Kong #endif 9808f3b4b4dSDmitry Karpeev if (iis) osm->user_subdomains = PETSC_TRUE; 981b1a0a8a3SJed Brown PetscFunctionReturn(0); 982b1a0a8a3SJed Brown } 983b1a0a8a3SJed Brown 984b1a0a8a3SJed Brown 985f7a08781SBarry Smith static PetscErrorCode PCGASMSetOverlap_GASM(PC pc,PetscInt ovl) 986b1a0a8a3SJed Brown { 987f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 988b1a0a8a3SJed Brown 989b1a0a8a3SJed Brown PetscFunctionBegin; 990ce94432eSBarry Smith if (ovl < 0) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap value requested"); 991ce94432eSBarry Smith if (pc->setupcalled && ovl != osm->overlap) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetOverlap() should be called before PCSetUp()."); 9922fa5cd67SKarl Rupp if (!pc->setupcalled) osm->overlap = ovl; 993b1a0a8a3SJed Brown PetscFunctionReturn(0); 994b1a0a8a3SJed Brown } 995b1a0a8a3SJed Brown 996f7a08781SBarry Smith static PetscErrorCode PCGASMSetType_GASM(PC pc,PCGASMType type) 997b1a0a8a3SJed Brown { 998f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 999b1a0a8a3SJed Brown 1000b1a0a8a3SJed Brown PetscFunctionBegin; 1001b1a0a8a3SJed Brown osm->type = type; 1002b1a0a8a3SJed Brown osm->type_set = PETSC_TRUE; 1003b1a0a8a3SJed Brown PetscFunctionReturn(0); 1004b1a0a8a3SJed Brown } 1005b1a0a8a3SJed Brown 1006f7a08781SBarry Smith static PetscErrorCode PCGASMSetSortIndices_GASM(PC pc,PetscBool doSort) 1007b1a0a8a3SJed Brown { 1008f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 1009b1a0a8a3SJed Brown 1010b1a0a8a3SJed Brown PetscFunctionBegin; 1011b1a0a8a3SJed Brown osm->sort_indices = doSort; 1012b1a0a8a3SJed Brown PetscFunctionReturn(0); 1013b1a0a8a3SJed Brown } 1014b1a0a8a3SJed Brown 101544fe18e5SDmitry Karpeev /* 10168f3b4b4dSDmitry Karpeev FIXME: This routine might need to be modified now that multiple ranks per subdomain are allowed. 101744fe18e5SDmitry Karpeev In particular, it would upset the global subdomain number calculation. 101844fe18e5SDmitry Karpeev */ 1019f7a08781SBarry Smith static PetscErrorCode PCGASMGetSubKSP_GASM(PC pc,PetscInt *n,PetscInt *first,KSP **ksp) 1020b1a0a8a3SJed Brown { 1021f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 1022b1a0a8a3SJed Brown PetscErrorCode ierr; 1023b1a0a8a3SJed Brown 1024b1a0a8a3SJed Brown PetscFunctionBegin; 1025ce94432eSBarry Smith if (osm->n < 1) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Need to call PCSetUP() on PC (or KSPSetUp() on the outer KSP object) before calling here"); 1026b1a0a8a3SJed Brown 10272fa5cd67SKarl Rupp if (n) *n = osm->n; 1028ab3e923fSDmitry Karpeev if (first) { 1029ce94432eSBarry Smith ierr = MPI_Scan(&osm->n,first,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr); 1030ab3e923fSDmitry Karpeev *first -= osm->n; 1031b1a0a8a3SJed Brown } 1032b1a0a8a3SJed Brown if (ksp) { 1033b1a0a8a3SJed Brown /* Assume that local solves are now different; not necessarily 103406b43e7eSDmitry Karpeev true, though! This flag is used only for PCView_GASM() */ 1035b1a0a8a3SJed Brown *ksp = osm->ksp; 10366a4f0f73SDmitry Karpeev osm->same_subdomain_solvers = PETSC_FALSE; 1037b1a0a8a3SJed Brown } 1038b1a0a8a3SJed Brown PetscFunctionReturn(0); 1039ab3e923fSDmitry Karpeev } /* PCGASMGetSubKSP_GASM() */ 1040b1a0a8a3SJed Brown 1041b1a0a8a3SJed Brown /*@C 104206b43e7eSDmitry Karpeev PCGASMSetSubdomains - Sets the subdomains for this processor 104306b43e7eSDmitry Karpeev for the additive Schwarz preconditioner. 1044b1a0a8a3SJed Brown 10458f3b4b4dSDmitry Karpeev Collective on pc 1046b1a0a8a3SJed Brown 1047b1a0a8a3SJed Brown Input Parameters: 10488f3b4b4dSDmitry Karpeev + pc - the preconditioner object 104906b43e7eSDmitry Karpeev . n - the number of subdomains for this processor 10508f3b4b4dSDmitry Karpeev . iis - the index sets that define the inner subdomains (or NULL for PETSc to determine subdomains) 10518f3b4b4dSDmitry 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) 1052b1a0a8a3SJed Brown 1053b1a0a8a3SJed Brown Notes: 105406b43e7eSDmitry Karpeev The IS indices use the parallel, global numbering of the vector entries. 10556a4f0f73SDmitry Karpeev Inner subdomains are those where the correction is applied. 10566a4f0f73SDmitry Karpeev Outer subdomains are those where the residual necessary to obtain the 10576a4f0f73SDmitry Karpeev corrections is obtained (see PCGASMType for the use of inner/outer subdomains). 10586a4f0f73SDmitry Karpeev Both inner and outer subdomains can extend over several processors. 10596a4f0f73SDmitry Karpeev This processor's portion of a subdomain is known as a local subdomain. 10606a4f0f73SDmitry Karpeev 1061670417b2SFande Kong Inner subdomains can not overlap with each other, do not have any entities from remote processors, 1062670417b2SFande Kong and have to cover the entire local subdomain owned by the current processor. The index sets on each 1063670417b2SFande Kong process should be ordered such that the ith local subdomain is connected to the ith remote subdomain 1064670417b2SFande Kong on another MPI process. 1065670417b2SFande Kong 10666a4f0f73SDmitry Karpeev By default the GASM preconditioner uses 1 (local) subdomain per processor. 10676a4f0f73SDmitry Karpeev 1068b1a0a8a3SJed Brown 1069b1a0a8a3SJed Brown Level: advanced 1070b1a0a8a3SJed Brown 107106b43e7eSDmitry Karpeev .keywords: PC, GASM, set, subdomains, additive Schwarz 1072b1a0a8a3SJed Brown 10738f3b4b4dSDmitry Karpeev .seealso: PCGASMSetNumSubdomains(), PCGASMSetOverlap(), PCGASMGetSubKSP(), 107406b43e7eSDmitry Karpeev PCGASMCreateSubdomains2D(), PCGASMGetSubdomains() 1075b1a0a8a3SJed Brown @*/ 10766a4f0f73SDmitry Karpeev PetscErrorCode PCGASMSetSubdomains(PC pc,PetscInt n,IS iis[],IS ois[]) 1077b1a0a8a3SJed Brown { 10788f3b4b4dSDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 1079b1a0a8a3SJed Brown PetscErrorCode ierr; 1080b1a0a8a3SJed Brown 1081b1a0a8a3SJed Brown PetscFunctionBegin; 1082b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 10836a4f0f73SDmitry Karpeev ierr = PetscTryMethod(pc,"PCGASMSetSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,n,iis,ois));CHKERRQ(ierr); 10848f3b4b4dSDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 1085b1a0a8a3SJed Brown PetscFunctionReturn(0); 1086b1a0a8a3SJed Brown } 1087b1a0a8a3SJed Brown 1088b1a0a8a3SJed Brown 1089b1a0a8a3SJed Brown /*@ 1090f746d493SDmitry Karpeev PCGASMSetOverlap - Sets the overlap between a pair of subdomains for the 1091b1a0a8a3SJed Brown additive Schwarz preconditioner. Either all or no processors in the 10928f3b4b4dSDmitry Karpeev pc communicator must call this routine. 1093b1a0a8a3SJed Brown 10948f3b4b4dSDmitry Karpeev Logically Collective on pc 1095b1a0a8a3SJed Brown 1096b1a0a8a3SJed Brown Input Parameters: 1097b1a0a8a3SJed Brown + pc - the preconditioner context 10988f3b4b4dSDmitry Karpeev - ovl - the amount of overlap between subdomains (ovl >= 0, default value = 0) 1099b1a0a8a3SJed Brown 1100b1a0a8a3SJed Brown Options Database Key: 110106b43e7eSDmitry Karpeev . -pc_gasm_overlap <overlap> - Sets overlap 1102b1a0a8a3SJed Brown 1103b1a0a8a3SJed Brown Notes: 110406b43e7eSDmitry Karpeev By default the GASM preconditioner uses 1 subdomain per processor. To use 11058f3b4b4dSDmitry Karpeev multiple subdomain per perocessor or "straddling" subdomains that intersect 11068f3b4b4dSDmitry Karpeev multiple processors use PCGASMSetSubdomains() (or option -pc_gasm_total_subdomains <n>). 1107b1a0a8a3SJed Brown 11088f3b4b4dSDmitry Karpeev The overlap defaults to 0, so if one desires that no additional 1109b1a0a8a3SJed Brown overlap be computed beyond what may have been set with a call to 11108f3b4b4dSDmitry Karpeev PCGASMSetSubdomains(), then ovl must be set to be 0. In particular, if one does 11118f3b4b4dSDmitry Karpeev not explicitly set the subdomains in application code, then all overlap would be computed 1112f746d493SDmitry Karpeev internally by PETSc, and using an overlap of 0 would result in an GASM 1113b1a0a8a3SJed Brown variant that is equivalent to the block Jacobi preconditioner. 1114b1a0a8a3SJed Brown 1115b1a0a8a3SJed Brown Note that one can define initial index sets with any overlap via 111606b43e7eSDmitry Karpeev PCGASMSetSubdomains(); the routine PCGASMSetOverlap() merely allows 111706b43e7eSDmitry Karpeev PETSc to extend that overlap further, if desired. 1118b1a0a8a3SJed Brown 1119b1a0a8a3SJed Brown Level: intermediate 1120b1a0a8a3SJed Brown 1121f746d493SDmitry Karpeev .keywords: PC, GASM, set, overlap 1122b1a0a8a3SJed Brown 11238f3b4b4dSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMGetSubKSP(), 112406b43e7eSDmitry Karpeev PCGASMCreateSubdomains2D(), PCGASMGetSubdomains() 1125b1a0a8a3SJed Brown @*/ 11267087cfbeSBarry Smith PetscErrorCode PCGASMSetOverlap(PC pc,PetscInt ovl) 1127b1a0a8a3SJed Brown { 1128b1a0a8a3SJed Brown PetscErrorCode ierr; 11298f3b4b4dSDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 1130b1a0a8a3SJed Brown 1131b1a0a8a3SJed Brown PetscFunctionBegin; 1132b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1133b1a0a8a3SJed Brown PetscValidLogicalCollectiveInt(pc,ovl,2); 1134f746d493SDmitry Karpeev ierr = PetscTryMethod(pc,"PCGASMSetOverlap_C",(PC,PetscInt),(pc,ovl));CHKERRQ(ierr); 11358f3b4b4dSDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 1136b1a0a8a3SJed Brown PetscFunctionReturn(0); 1137b1a0a8a3SJed Brown } 1138b1a0a8a3SJed Brown 1139b1a0a8a3SJed Brown /*@ 1140f746d493SDmitry Karpeev PCGASMSetType - Sets the type of restriction and interpolation used 1141b1a0a8a3SJed Brown for local problems in the additive Schwarz method. 1142b1a0a8a3SJed Brown 1143b1a0a8a3SJed Brown Logically Collective on PC 1144b1a0a8a3SJed Brown 1145b1a0a8a3SJed Brown Input Parameters: 1146b1a0a8a3SJed Brown + pc - the preconditioner context 1147f746d493SDmitry Karpeev - type - variant of GASM, one of 1148b1a0a8a3SJed Brown .vb 1149f746d493SDmitry Karpeev PC_GASM_BASIC - full interpolation and restriction 1150f746d493SDmitry Karpeev PC_GASM_RESTRICT - full restriction, local processor interpolation 1151f746d493SDmitry Karpeev PC_GASM_INTERPOLATE - full interpolation, local processor restriction 1152f746d493SDmitry Karpeev PC_GASM_NONE - local processor restriction and interpolation 1153b1a0a8a3SJed Brown .ve 1154b1a0a8a3SJed Brown 1155b1a0a8a3SJed Brown Options Database Key: 1156f746d493SDmitry Karpeev . -pc_gasm_type [basic,restrict,interpolate,none] - Sets GASM type 1157b1a0a8a3SJed Brown 1158b1a0a8a3SJed Brown Level: intermediate 1159b1a0a8a3SJed Brown 1160f746d493SDmitry Karpeev .keywords: PC, GASM, set, type 1161b1a0a8a3SJed Brown 11628f3b4b4dSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMGetSubKSP(), 1163f746d493SDmitry Karpeev PCGASMCreateSubdomains2D() 1164b1a0a8a3SJed Brown @*/ 11657087cfbeSBarry Smith PetscErrorCode PCGASMSetType(PC pc,PCGASMType type) 1166b1a0a8a3SJed Brown { 1167b1a0a8a3SJed Brown PetscErrorCode ierr; 1168b1a0a8a3SJed Brown 1169b1a0a8a3SJed Brown PetscFunctionBegin; 1170b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1171b1a0a8a3SJed Brown PetscValidLogicalCollectiveEnum(pc,type,2); 1172f746d493SDmitry Karpeev ierr = PetscTryMethod(pc,"PCGASMSetType_C",(PC,PCGASMType),(pc,type));CHKERRQ(ierr); 1173b1a0a8a3SJed Brown PetscFunctionReturn(0); 1174b1a0a8a3SJed Brown } 1175b1a0a8a3SJed Brown 1176b1a0a8a3SJed Brown /*@ 1177f746d493SDmitry Karpeev PCGASMSetSortIndices - Determines whether subdomain indices are sorted. 1178b1a0a8a3SJed Brown 1179b1a0a8a3SJed Brown Logically Collective on PC 1180b1a0a8a3SJed Brown 1181b1a0a8a3SJed Brown Input Parameters: 1182b1a0a8a3SJed Brown + pc - the preconditioner context 1183b1a0a8a3SJed Brown - doSort - sort the subdomain indices 1184b1a0a8a3SJed Brown 1185b1a0a8a3SJed Brown Level: intermediate 1186b1a0a8a3SJed Brown 1187f746d493SDmitry Karpeev .keywords: PC, GASM, set, type 1188b1a0a8a3SJed Brown 11898f3b4b4dSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMGetSubKSP(), 1190f746d493SDmitry Karpeev PCGASMCreateSubdomains2D() 1191b1a0a8a3SJed Brown @*/ 11927087cfbeSBarry Smith PetscErrorCode PCGASMSetSortIndices(PC pc,PetscBool doSort) 1193b1a0a8a3SJed Brown { 1194b1a0a8a3SJed Brown PetscErrorCode ierr; 1195b1a0a8a3SJed Brown 1196b1a0a8a3SJed Brown PetscFunctionBegin; 1197b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1198b1a0a8a3SJed Brown PetscValidLogicalCollectiveBool(pc,doSort,2); 1199f746d493SDmitry Karpeev ierr = PetscTryMethod(pc,"PCGASMSetSortIndices_C",(PC,PetscBool),(pc,doSort));CHKERRQ(ierr); 1200b1a0a8a3SJed Brown PetscFunctionReturn(0); 1201b1a0a8a3SJed Brown } 1202b1a0a8a3SJed Brown 1203b1a0a8a3SJed Brown /*@C 1204f746d493SDmitry Karpeev PCGASMGetSubKSP - Gets the local KSP contexts for all blocks on 1205b1a0a8a3SJed Brown this processor. 1206b1a0a8a3SJed Brown 1207b1a0a8a3SJed Brown Collective on PC iff first_local is requested 1208b1a0a8a3SJed Brown 1209b1a0a8a3SJed Brown Input Parameter: 1210b1a0a8a3SJed Brown . pc - the preconditioner context 1211b1a0a8a3SJed Brown 1212b1a0a8a3SJed Brown Output Parameters: 12130298fd71SBarry Smith + n_local - the number of blocks on this processor or NULL 12140298fd71SBarry Smith . first_local - the global number of the first block on this processor or NULL, 12150298fd71SBarry Smith all processors must request or all must pass NULL 1216b1a0a8a3SJed Brown - ksp - the array of KSP contexts 1217b1a0a8a3SJed Brown 1218b1a0a8a3SJed Brown Note: 1219f746d493SDmitry Karpeev After PCGASMGetSubKSP() the array of KSPes is not to be freed 1220b1a0a8a3SJed Brown 1221b1a0a8a3SJed Brown Currently for some matrix implementations only 1 block per processor 1222b1a0a8a3SJed Brown is supported. 1223b1a0a8a3SJed Brown 1224f746d493SDmitry Karpeev You must call KSPSetUp() before calling PCGASMGetSubKSP(). 1225b1a0a8a3SJed Brown 1226b1a0a8a3SJed Brown Level: advanced 1227b1a0a8a3SJed Brown 1228f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, get, sub, KSP, context 1229b1a0a8a3SJed Brown 12308f3b4b4dSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMSetOverlap(), 1231f746d493SDmitry Karpeev PCGASMCreateSubdomains2D(), 1232b1a0a8a3SJed Brown @*/ 12337087cfbeSBarry Smith PetscErrorCode PCGASMGetSubKSP(PC pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp[]) 1234b1a0a8a3SJed Brown { 1235b1a0a8a3SJed Brown PetscErrorCode ierr; 1236b1a0a8a3SJed Brown 1237b1a0a8a3SJed Brown PetscFunctionBegin; 1238b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1239f746d493SDmitry Karpeev ierr = PetscUseMethod(pc,"PCGASMGetSubKSP_C",(PC,PetscInt*,PetscInt*,KSP **),(pc,n_local,first_local,ksp));CHKERRQ(ierr); 1240b1a0a8a3SJed Brown PetscFunctionReturn(0); 1241b1a0a8a3SJed Brown } 1242b1a0a8a3SJed Brown 1243b1a0a8a3SJed Brown /* -------------------------------------------------------------------------------------*/ 1244b1a0a8a3SJed Brown /*MC 1245f746d493SDmitry Karpeev PCGASM - Use the (restricted) additive Schwarz method, each block is (approximately) solved with 1246b1a0a8a3SJed Brown its own KSP object. 1247b1a0a8a3SJed Brown 1248b1a0a8a3SJed Brown Options Database Keys: 12498f3b4b4dSDmitry Karpeev + -pc_gasm_total_subdomains <n> - Sets total number of local subdomains to be distributed among processors 125006b43e7eSDmitry Karpeev . -pc_gasm_view_subdomains - activates the printing of subdomain indices in PCView(), -ksp_view or -snes_view 125106b43e7eSDmitry Karpeev . -pc_gasm_print_subdomains - activates the printing of subdomain indices in PCSetUp() 125206b43e7eSDmitry Karpeev . -pc_gasm_overlap <ovl> - Sets overlap by which to (automatically) extend local subdomains 1253f746d493SDmitry Karpeev - -pc_gasm_type [basic,restrict,interpolate,none] - Sets GASM type 1254b1a0a8a3SJed Brown 1255b1a0a8a3SJed Brown IMPORTANT: If you run with, for example, 3 blocks on 1 processor or 3 blocks on 3 processors you 1256f746d493SDmitry Karpeev will get a different convergence rate due to the default option of -pc_gasm_type restrict. Use 1257f746d493SDmitry Karpeev -pc_gasm_type basic to use the standard GASM. 1258b1a0a8a3SJed Brown 125995452b02SPatrick Sanan Notes: 126095452b02SPatrick Sanan Blocks can be shared by multiple processes. 1261b1a0a8a3SJed Brown 1262b1a0a8a3SJed Brown To set options on the solvers for each block append -sub_ to all the KSP, and PC 1263b1a0a8a3SJed Brown options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly 1264b1a0a8a3SJed Brown 1265f746d493SDmitry Karpeev To set the options on the solvers separate for each block call PCGASMGetSubKSP() 1266b1a0a8a3SJed Brown and set the options directly on the resulting KSP object (you can access its PC 1267b1a0a8a3SJed Brown with KSPGetPC()) 1268b1a0a8a3SJed Brown 1269b1a0a8a3SJed Brown 1270b1a0a8a3SJed Brown Level: beginner 1271b1a0a8a3SJed Brown 1272b1a0a8a3SJed Brown Concepts: additive Schwarz method 1273b1a0a8a3SJed Brown 1274b1a0a8a3SJed Brown References: 127596a0c994SBarry Smith + 1. - M Dryja, OB Widlund, An additive variant of the Schwarz alternating method for the case of many subregions 127696a0c994SBarry Smith Courant Institute, New York University Technical report 127796a0c994SBarry Smith - 2. - Barry Smith, Petter Bjorstad, and William Gropp, Domain Decompositions: Parallel Multilevel Methods for Elliptic Partial Differential Equations, 127896a0c994SBarry Smith Cambridge University Press. 1279b1a0a8a3SJed Brown 1280b1a0a8a3SJed Brown .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 128149517cdeSBarry Smith PCBJACOBI, PCGASMGetSubKSP(), PCGASMSetSubdomains(), 12828f3b4b4dSDmitry Karpeev PCSetModifySubmatrices(), PCGASMSetOverlap(), PCGASMSetType() 1283b1a0a8a3SJed Brown 1284b1a0a8a3SJed Brown M*/ 1285b1a0a8a3SJed Brown 12868cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_GASM(PC pc) 1287b1a0a8a3SJed Brown { 1288b1a0a8a3SJed Brown PetscErrorCode ierr; 1289f746d493SDmitry Karpeev PC_GASM *osm; 1290b1a0a8a3SJed Brown 1291b1a0a8a3SJed Brown PetscFunctionBegin; 1292b00a9115SJed Brown ierr = PetscNewLog(pc,&osm);CHKERRQ(ierr); 12932fa5cd67SKarl Rupp 12948f3b4b4dSDmitry Karpeev osm->N = PETSC_DETERMINE; 129506b43e7eSDmitry Karpeev osm->n = PETSC_DECIDE; 12968f3b4b4dSDmitry Karpeev osm->nmax = PETSC_DETERMINE; 12978f3b4b4dSDmitry Karpeev osm->overlap = 0; 1298b1a0a8a3SJed Brown osm->ksp = 0; 12996a4f0f73SDmitry Karpeev osm->gorestriction = 0; 13006a4f0f73SDmitry Karpeev osm->girestriction = 0; 1301ea91fabdSFande Kong osm->pctoouter = 0; 1302ab3e923fSDmitry Karpeev osm->gx = 0; 1303ab3e923fSDmitry Karpeev osm->gy = 0; 1304b1a0a8a3SJed Brown osm->x = 0; 1305b1a0a8a3SJed Brown osm->y = 0; 1306ea91fabdSFande Kong osm->pcx = 0; 1307ea91fabdSFande Kong osm->pcy = 0; 1308ea91fabdSFande Kong osm->permutationIS = 0; 1309ea91fabdSFande Kong osm->permutationP = 0; 1310ea91fabdSFande Kong osm->pcmat = 0; 13116a4f0f73SDmitry Karpeev osm->ois = 0; 13126a4f0f73SDmitry Karpeev osm->iis = 0; 1313b1a0a8a3SJed Brown osm->pmat = 0; 1314f746d493SDmitry Karpeev osm->type = PC_GASM_RESTRICT; 13156a4f0f73SDmitry Karpeev osm->same_subdomain_solvers = PETSC_TRUE; 1316b1a0a8a3SJed Brown osm->sort_indices = PETSC_TRUE; 1317d709ea83SDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 1318ea91fabdSFande Kong osm->hierarchicalpartitioning = PETSC_FALSE; 1319b1a0a8a3SJed Brown 1320b1a0a8a3SJed Brown pc->data = (void*)osm; 1321f746d493SDmitry Karpeev pc->ops->apply = PCApply_GASM; 1322f746d493SDmitry Karpeev pc->ops->applytranspose = PCApplyTranspose_GASM; 1323f746d493SDmitry Karpeev pc->ops->setup = PCSetUp_GASM; 1324a06653b4SBarry Smith pc->ops->reset = PCReset_GASM; 1325f746d493SDmitry Karpeev pc->ops->destroy = PCDestroy_GASM; 1326f746d493SDmitry Karpeev pc->ops->setfromoptions = PCSetFromOptions_GASM; 1327f746d493SDmitry Karpeev pc->ops->setuponblocks = PCSetUpOnBlocks_GASM; 1328f746d493SDmitry Karpeev pc->ops->view = PCView_GASM; 1329b1a0a8a3SJed Brown pc->ops->applyrichardson = 0; 1330b1a0a8a3SJed Brown 1331bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetSubdomains_C",PCGASMSetSubdomains_GASM);CHKERRQ(ierr); 1332bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetOverlap_C",PCGASMSetOverlap_GASM);CHKERRQ(ierr); 1333bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetType_C",PCGASMSetType_GASM);CHKERRQ(ierr); 1334bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetSortIndices_C",PCGASMSetSortIndices_GASM);CHKERRQ(ierr); 1335bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGASMGetSubKSP_C",PCGASMGetSubKSP_GASM);CHKERRQ(ierr); 1336b1a0a8a3SJed Brown PetscFunctionReturn(0); 1337b1a0a8a3SJed Brown } 1338b1a0a8a3SJed Brown 1339b1a0a8a3SJed Brown 13408f3b4b4dSDmitry Karpeev PetscErrorCode PCGASMCreateLocalSubdomains(Mat A, PetscInt nloc, IS *iis[]) 1341b1a0a8a3SJed Brown { 1342b1a0a8a3SJed Brown MatPartitioning mpart; 1343b1a0a8a3SJed Brown const char *prefix; 1344b1a0a8a3SJed Brown PetscInt i,j,rstart,rend,bs; 1345976e8c5aSLisandro Dalcin PetscBool hasop, isbaij = PETSC_FALSE,foundpart = PETSC_FALSE; 13460298fd71SBarry Smith Mat Ad = NULL, adj; 1347b1a0a8a3SJed Brown IS ispart,isnumb,*is; 1348b1a0a8a3SJed Brown PetscErrorCode ierr; 1349b1a0a8a3SJed Brown 1350b1a0a8a3SJed Brown PetscFunctionBegin; 13518f3b4b4dSDmitry Karpeev if (nloc < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of local subdomains must > 0, got nloc = %D",nloc); 1352b1a0a8a3SJed Brown 1353b1a0a8a3SJed Brown /* Get prefix, row distribution, and block size */ 1354b1a0a8a3SJed Brown ierr = MatGetOptionsPrefix(A,&prefix);CHKERRQ(ierr); 1355b1a0a8a3SJed Brown ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 1356b1a0a8a3SJed Brown ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr); 1357b1a0a8a3SJed 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); 1358b1a0a8a3SJed Brown 1359b1a0a8a3SJed Brown /* Get diagonal block from matrix if possible */ 1360976e8c5aSLisandro Dalcin ierr = MatHasOperation(A,MATOP_GET_DIAGONAL_BLOCK,&hasop);CHKERRQ(ierr); 1361976e8c5aSLisandro Dalcin if (hasop) { 136211bd1e4dSLisandro Dalcin ierr = MatGetDiagonalBlock(A,&Ad);CHKERRQ(ierr); 1363b1a0a8a3SJed Brown } 1364b1a0a8a3SJed Brown if (Ad) { 1365251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)Ad,MATSEQBAIJ,&isbaij);CHKERRQ(ierr); 1366251f4c67SDmitry Karpeev if (!isbaij) {ierr = PetscObjectTypeCompare((PetscObject)Ad,MATSEQSBAIJ,&isbaij);CHKERRQ(ierr);} 1367b1a0a8a3SJed Brown } 13688f3b4b4dSDmitry Karpeev if (Ad && nloc > 1) { 1369b1a0a8a3SJed Brown PetscBool match,done; 1370b1a0a8a3SJed Brown /* Try to setup a good matrix partitioning if available */ 1371b1a0a8a3SJed Brown ierr = MatPartitioningCreate(PETSC_COMM_SELF,&mpart);CHKERRQ(ierr); 1372b1a0a8a3SJed Brown ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr); 1373b1a0a8a3SJed Brown ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr); 1374251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGCURRENT,&match);CHKERRQ(ierr); 1375b1a0a8a3SJed Brown if (!match) { 1376251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGSQUARE,&match);CHKERRQ(ierr); 1377b1a0a8a3SJed Brown } 1378b1a0a8a3SJed Brown if (!match) { /* assume a "good" partitioner is available */ 13791a83f524SJed Brown PetscInt na; 13801a83f524SJed Brown const PetscInt *ia,*ja; 1381b1a0a8a3SJed Brown ierr = MatGetRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr); 1382b1a0a8a3SJed Brown if (done) { 1383b1a0a8a3SJed Brown /* Build adjacency matrix by hand. Unfortunately a call to 1384b1a0a8a3SJed Brown MatConvert(Ad,MATMPIADJ,MAT_INITIAL_MATRIX,&adj) will 1385b1a0a8a3SJed Brown remove the block-aij structure and we cannot expect 1386b1a0a8a3SJed Brown MatPartitioning to split vertices as we need */ 13871a83f524SJed Brown PetscInt i,j,len,nnz,cnt,*iia=0,*jja=0; 13881a83f524SJed Brown const PetscInt *row; 1389b1a0a8a3SJed Brown nnz = 0; 1390b1a0a8a3SJed Brown for (i=0; i<na; i++) { /* count number of nonzeros */ 1391b1a0a8a3SJed Brown len = ia[i+1] - ia[i]; 1392b1a0a8a3SJed Brown row = ja + ia[i]; 1393b1a0a8a3SJed Brown for (j=0; j<len; j++) { 1394b1a0a8a3SJed Brown if (row[j] == i) { /* don't count diagonal */ 1395b1a0a8a3SJed Brown len--; break; 1396b1a0a8a3SJed Brown } 1397b1a0a8a3SJed Brown } 1398b1a0a8a3SJed Brown nnz += len; 1399b1a0a8a3SJed Brown } 1400854ce69bSBarry Smith ierr = PetscMalloc1(na+1,&iia);CHKERRQ(ierr); 1401854ce69bSBarry Smith ierr = PetscMalloc1(nnz,&jja);CHKERRQ(ierr); 1402b1a0a8a3SJed Brown nnz = 0; 1403b1a0a8a3SJed Brown iia[0] = 0; 1404b1a0a8a3SJed Brown for (i=0; i<na; i++) { /* fill adjacency */ 1405b1a0a8a3SJed Brown cnt = 0; 1406b1a0a8a3SJed Brown len = ia[i+1] - ia[i]; 1407b1a0a8a3SJed Brown row = ja + ia[i]; 1408b1a0a8a3SJed Brown for (j=0; j<len; j++) { 14092fa5cd67SKarl Rupp if (row[j] != i) jja[nnz+cnt++] = row[j]; /* if not diagonal */ 1410b1a0a8a3SJed Brown } 1411b1a0a8a3SJed Brown nnz += cnt; 1412b1a0a8a3SJed Brown iia[i+1] = nnz; 1413b1a0a8a3SJed Brown } 1414b1a0a8a3SJed Brown /* Partitioning of the adjacency matrix */ 14150298fd71SBarry Smith ierr = MatCreateMPIAdj(PETSC_COMM_SELF,na,na,iia,jja,NULL,&adj);CHKERRQ(ierr); 1416b1a0a8a3SJed Brown ierr = MatPartitioningSetAdjacency(mpart,adj);CHKERRQ(ierr); 14178f3b4b4dSDmitry Karpeev ierr = MatPartitioningSetNParts(mpart,nloc);CHKERRQ(ierr); 1418b1a0a8a3SJed Brown ierr = MatPartitioningApply(mpart,&ispart);CHKERRQ(ierr); 1419b1a0a8a3SJed Brown ierr = ISPartitioningToNumbering(ispart,&isnumb);CHKERRQ(ierr); 14206bf464f9SBarry Smith ierr = MatDestroy(&adj);CHKERRQ(ierr); 1421b1a0a8a3SJed Brown foundpart = PETSC_TRUE; 1422b1a0a8a3SJed Brown } 1423b1a0a8a3SJed Brown ierr = MatRestoreRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr); 1424b1a0a8a3SJed Brown } 14256bf464f9SBarry Smith ierr = MatPartitioningDestroy(&mpart);CHKERRQ(ierr); 1426b1a0a8a3SJed Brown } 14278f3b4b4dSDmitry Karpeev ierr = PetscMalloc1(nloc,&is);CHKERRQ(ierr); 1428b1a0a8a3SJed Brown if (!foundpart) { 1429b1a0a8a3SJed Brown 1430b1a0a8a3SJed Brown /* Partitioning by contiguous chunks of rows */ 1431b1a0a8a3SJed Brown 1432b1a0a8a3SJed Brown PetscInt mbs = (rend-rstart)/bs; 1433b1a0a8a3SJed Brown PetscInt start = rstart; 14348f3b4b4dSDmitry Karpeev for (i=0; i<nloc; i++) { 14358f3b4b4dSDmitry Karpeev PetscInt count = (mbs/nloc + ((mbs % nloc) > i)) * bs; 1436b1a0a8a3SJed Brown ierr = ISCreateStride(PETSC_COMM_SELF,count,start,1,&is[i]);CHKERRQ(ierr); 1437b1a0a8a3SJed Brown start += count; 1438b1a0a8a3SJed Brown } 1439b1a0a8a3SJed Brown 1440b1a0a8a3SJed Brown } else { 1441b1a0a8a3SJed Brown 1442b1a0a8a3SJed Brown /* Partitioning by adjacency of diagonal block */ 1443b1a0a8a3SJed Brown 1444b1a0a8a3SJed Brown const PetscInt *numbering; 1445b1a0a8a3SJed Brown PetscInt *count,nidx,*indices,*newidx,start=0; 1446b1a0a8a3SJed Brown /* Get node count in each partition */ 14478f3b4b4dSDmitry Karpeev ierr = PetscMalloc1(nloc,&count);CHKERRQ(ierr); 14488f3b4b4dSDmitry Karpeev ierr = ISPartitioningCount(ispart,nloc,count);CHKERRQ(ierr); 1449b1a0a8a3SJed Brown if (isbaij && bs > 1) { /* adjust for the block-aij case */ 14508f3b4b4dSDmitry Karpeev for (i=0; i<nloc; i++) count[i] *= bs; 1451b1a0a8a3SJed Brown } 1452b1a0a8a3SJed Brown /* Build indices from node numbering */ 1453b1a0a8a3SJed Brown ierr = ISGetLocalSize(isnumb,&nidx);CHKERRQ(ierr); 1454785e854fSJed Brown ierr = PetscMalloc1(nidx,&indices);CHKERRQ(ierr); 1455b1a0a8a3SJed Brown for (i=0; i<nidx; i++) indices[i] = i; /* needs to be initialized */ 1456b1a0a8a3SJed Brown ierr = ISGetIndices(isnumb,&numbering);CHKERRQ(ierr); 1457b1a0a8a3SJed Brown ierr = PetscSortIntWithPermutation(nidx,numbering,indices);CHKERRQ(ierr); 1458b1a0a8a3SJed Brown ierr = ISRestoreIndices(isnumb,&numbering);CHKERRQ(ierr); 1459b1a0a8a3SJed Brown if (isbaij && bs > 1) { /* adjust for the block-aij case */ 1460785e854fSJed Brown ierr = PetscMalloc1(nidx*bs,&newidx);CHKERRQ(ierr); 14612fa5cd67SKarl Rupp for (i=0; i<nidx; i++) { 14622fa5cd67SKarl Rupp for (j=0; j<bs; j++) newidx[i*bs+j] = indices[i]*bs + j; 14632fa5cd67SKarl Rupp } 1464b1a0a8a3SJed Brown ierr = PetscFree(indices);CHKERRQ(ierr); 1465b1a0a8a3SJed Brown nidx *= bs; 1466b1a0a8a3SJed Brown indices = newidx; 1467b1a0a8a3SJed Brown } 1468b1a0a8a3SJed Brown /* Shift to get global indices */ 1469b1a0a8a3SJed Brown for (i=0; i<nidx; i++) indices[i] += rstart; 1470b1a0a8a3SJed Brown 1471b1a0a8a3SJed Brown /* Build the index sets for each block */ 14728f3b4b4dSDmitry Karpeev for (i=0; i<nloc; i++) { 1473b1a0a8a3SJed Brown ierr = ISCreateGeneral(PETSC_COMM_SELF,count[i],&indices[start],PETSC_COPY_VALUES,&is[i]);CHKERRQ(ierr); 1474b1a0a8a3SJed Brown ierr = ISSort(is[i]);CHKERRQ(ierr); 1475b1a0a8a3SJed Brown start += count[i]; 1476b1a0a8a3SJed Brown } 1477b1a0a8a3SJed Brown 1478302440fdSBarry Smith ierr = PetscFree(count);CHKERRQ(ierr); 1479302440fdSBarry Smith ierr = PetscFree(indices);CHKERRQ(ierr); 1480fcfd50ebSBarry Smith ierr = ISDestroy(&isnumb);CHKERRQ(ierr); 1481fcfd50ebSBarry Smith ierr = ISDestroy(&ispart);CHKERRQ(ierr); 1482b1a0a8a3SJed Brown } 14836a4f0f73SDmitry Karpeev *iis = is; 14848f3b4b4dSDmitry Karpeev PetscFunctionReturn(0); 14858f3b4b4dSDmitry Karpeev } 14868f3b4b4dSDmitry Karpeev 1487b541e6a4SDmitry Karpeev PETSC_INTERN PetscErrorCode PCGASMCreateStraddlingSubdomains(Mat A,PetscInt N,PetscInt *n,IS *iis[]) 14888f3b4b4dSDmitry Karpeev { 14898f3b4b4dSDmitry Karpeev PetscErrorCode ierr; 14908f3b4b4dSDmitry Karpeev 14918f3b4b4dSDmitry Karpeev PetscFunctionBegin; 14928f3b4b4dSDmitry Karpeev ierr = MatSubdomainsCreateCoalesce(A,N,n,iis);CHKERRQ(ierr); 14938f3b4b4dSDmitry Karpeev PetscFunctionReturn(0); 14948f3b4b4dSDmitry Karpeev } 14958f3b4b4dSDmitry Karpeev 14968f3b4b4dSDmitry Karpeev 14978f3b4b4dSDmitry Karpeev 14988f3b4b4dSDmitry Karpeev /*@C 14998f3b4b4dSDmitry Karpeev PCGASMCreateSubdomains - Creates n index sets defining n nonoverlapping subdomains for the additive 15008f3b4b4dSDmitry Karpeev Schwarz preconditioner for a any problem based on its matrix. 15018f3b4b4dSDmitry Karpeev 15028f3b4b4dSDmitry Karpeev Collective 15038f3b4b4dSDmitry Karpeev 15048f3b4b4dSDmitry Karpeev Input Parameters: 15058f3b4b4dSDmitry Karpeev + A - The global matrix operator 15068f3b4b4dSDmitry Karpeev - N - the number of global subdomains requested 15078f3b4b4dSDmitry Karpeev 15088f3b4b4dSDmitry Karpeev Output Parameters: 15098f3b4b4dSDmitry Karpeev + n - the number of subdomains created on this processor 15108f3b4b4dSDmitry Karpeev - iis - the array of index sets defining the local inner subdomains (on which the correction is applied) 15118f3b4b4dSDmitry Karpeev 15128f3b4b4dSDmitry Karpeev Level: advanced 15138f3b4b4dSDmitry Karpeev 15148f3b4b4dSDmitry Karpeev Note: When N >= A's communicator size, each subdomain is local -- contained within a single processor. 15158f3b4b4dSDmitry Karpeev When N < size, the subdomains are 'straddling' (processor boundaries) and are no longer local. 15168f3b4b4dSDmitry Karpeev The resulting subdomains can be use in PCGASMSetSubdomains(pc,n,iss,NULL). The overlapping 15178f3b4b4dSDmitry Karpeev outer subdomains will be automatically generated from these according to the requested amount of 15188f3b4b4dSDmitry Karpeev overlap; this is currently supported only with local subdomains. 15198f3b4b4dSDmitry Karpeev 15208f3b4b4dSDmitry Karpeev 15218f3b4b4dSDmitry Karpeev .keywords: PC, GASM, additive Schwarz, create, subdomains, unstructured grid 15228f3b4b4dSDmitry Karpeev 15238f3b4b4dSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMDestroySubdomains() 15248f3b4b4dSDmitry Karpeev @*/ 1525b541e6a4SDmitry Karpeev PetscErrorCode PCGASMCreateSubdomains(Mat A,PetscInt N,PetscInt *n,IS *iis[]) 15268f3b4b4dSDmitry Karpeev { 15278f3b4b4dSDmitry Karpeev PetscMPIInt size; 15288f3b4b4dSDmitry Karpeev PetscErrorCode ierr; 15298f3b4b4dSDmitry Karpeev 15308f3b4b4dSDmitry Karpeev PetscFunctionBegin; 15318f3b4b4dSDmitry Karpeev PetscValidHeaderSpecific(A,MAT_CLASSID,1); 15328f3b4b4dSDmitry Karpeev PetscValidPointer(iis,4); 15338f3b4b4dSDmitry Karpeev 15348f3b4b4dSDmitry Karpeev if (N < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Number of subdomains must be > 0, N = %D",N); 15358f3b4b4dSDmitry Karpeev ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr); 15368f3b4b4dSDmitry Karpeev if (N >= size) { 15378f3b4b4dSDmitry Karpeev *n = N/size + (N%size); 15388f3b4b4dSDmitry Karpeev ierr = PCGASMCreateLocalSubdomains(A,*n,iis);CHKERRQ(ierr); 15396a4f0f73SDmitry Karpeev } else { 15408f3b4b4dSDmitry Karpeev ierr = PCGASMCreateStraddlingSubdomains(A,N,n,iis);CHKERRQ(ierr); 15416a4f0f73SDmitry Karpeev } 1542b1a0a8a3SJed Brown PetscFunctionReturn(0); 1543b1a0a8a3SJed Brown } 1544b1a0a8a3SJed Brown 1545b1a0a8a3SJed Brown /*@C 1546f746d493SDmitry Karpeev PCGASMDestroySubdomains - Destroys the index sets created with 15478f3b4b4dSDmitry Karpeev PCGASMCreateSubdomains() or PCGASMCreateSubdomains2D. Should be 154806b43e7eSDmitry Karpeev called after setting subdomains with PCGASMSetSubdomains(). 1549b1a0a8a3SJed Brown 1550b1a0a8a3SJed Brown Collective 1551b1a0a8a3SJed Brown 1552b1a0a8a3SJed Brown Input Parameters: 1553b1a0a8a3SJed Brown + n - the number of index sets 15546a4f0f73SDmitry Karpeev . iis - the array of inner subdomains, 15550298fd71SBarry Smith - ois - the array of outer subdomains, can be NULL 1556b1a0a8a3SJed Brown 15576a4f0f73SDmitry Karpeev Level: intermediate 15586a4f0f73SDmitry Karpeev 155995452b02SPatrick Sanan Notes: 156095452b02SPatrick Sanan this is merely a convenience subroutine that walks each list, 15612c112581SDmitry Karpeev destroys each IS on the list, and then frees the list. At the end the 15622c112581SDmitry Karpeev list pointers are set to NULL. 1563b1a0a8a3SJed Brown 1564f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, create, subdomains, unstructured grid 1565b1a0a8a3SJed Brown 15668f3b4b4dSDmitry Karpeev .seealso: PCGASMCreateSubdomains(), PCGASMSetSubdomains() 1567b1a0a8a3SJed Brown @*/ 15682c112581SDmitry Karpeev PetscErrorCode PCGASMDestroySubdomains(PetscInt n,IS **iis,IS **ois) 1569b1a0a8a3SJed Brown { 1570b1a0a8a3SJed Brown PetscInt i; 1571b1a0a8a3SJed Brown PetscErrorCode ierr; 15725fd66863SKarl Rupp 1573b1a0a8a3SJed Brown PetscFunctionBegin; 1574a35b7b57SDmitry Karpeev if (n <= 0) PetscFunctionReturn(0); 15756a4f0f73SDmitry Karpeev if (ois) { 15762c112581SDmitry Karpeev PetscValidPointer(ois,3); 15772c112581SDmitry Karpeev if (*ois) { 15782c112581SDmitry Karpeev PetscValidPointer(*ois,3); 1579a35b7b57SDmitry Karpeev for (i=0; i<n; i++) { 15802c112581SDmitry Karpeev ierr = ISDestroy(&(*ois)[i]);CHKERRQ(ierr); 1581a35b7b57SDmitry Karpeev } 1582135757f6SDmitry Karpeev ierr = PetscFree((*ois));CHKERRQ(ierr); 15832c112581SDmitry Karpeev } 1584b1a0a8a3SJed Brown } 1585b9d0fdaaSFande Kong if (iis) { 1586b9d0fdaaSFande Kong PetscValidPointer(iis,2); 1587b9d0fdaaSFande Kong if (*iis) { 1588b9d0fdaaSFande Kong PetscValidPointer(*iis,2); 1589b9d0fdaaSFande Kong for (i=0; i<n; i++) { 1590b9d0fdaaSFande Kong ierr = ISDestroy(&(*iis)[i]);CHKERRQ(ierr); 1591b9d0fdaaSFande Kong } 1592b9d0fdaaSFande Kong ierr = PetscFree((*iis));CHKERRQ(ierr); 1593b9d0fdaaSFande Kong } 1594b9d0fdaaSFande Kong } 1595b1a0a8a3SJed Brown PetscFunctionReturn(0); 1596b1a0a8a3SJed Brown } 1597b1a0a8a3SJed Brown 1598af538404SDmitry Karpeev 1599af538404SDmitry Karpeev #define PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,xleft_loc,ylow_loc,xright_loc,yhigh_loc,n) \ 1600af538404SDmitry Karpeev { \ 1601af538404SDmitry Karpeev PetscInt first_row = first/M, last_row = last/M+1; \ 1602af538404SDmitry Karpeev /* \ 1603af538404SDmitry Karpeev Compute ylow_loc and yhigh_loc so that (ylow_loc,xleft) and (yhigh_loc,xright) are the corners \ 1604af538404SDmitry Karpeev of the bounding box of the intersection of the subdomain with the local ownership range (local \ 1605af538404SDmitry Karpeev subdomain). \ 1606af538404SDmitry Karpeev Also compute xleft_loc and xright_loc as the lower and upper bounds on the first and last rows \ 1607af538404SDmitry Karpeev of the intersection. \ 1608af538404SDmitry Karpeev */ \ 1609af538404SDmitry Karpeev /* ylow_loc is the grid row containing the first element of the local sumbdomain */ \ 1610eec7e3faSDmitry Karpeev *ylow_loc = PetscMax(first_row,ylow); \ 1611af538404SDmitry Karpeev /* xleft_loc is the offset of first element of the local subdomain within its grid row (might actually be outside the local subdomain) */ \ 1612eec7e3faSDmitry Karpeev *xleft_loc = *ylow_loc==first_row ? PetscMax(first%M,xleft) : xleft; \ 1613af538404SDmitry Karpeev /* yhigh_loc is the grid row above the last local subdomain element */ \ 1614eec7e3faSDmitry Karpeev *yhigh_loc = PetscMin(last_row,yhigh); \ 1615af538404SDmitry Karpeev /* xright is the offset of the end of the local subdomain within its grid row (might actually be outside the local subdomain) */ \ 1616eec7e3faSDmitry Karpeev *xright_loc = *yhigh_loc==last_row ? PetscMin(xright,last%M) : xright; \ 1617af538404SDmitry Karpeev /* Now compute the size of the local subdomain n. */ \ 1618c3518bceSDmitry Karpeev *n = 0; \ 1619eec7e3faSDmitry Karpeev if (*ylow_loc < *yhigh_loc) { \ 1620af538404SDmitry Karpeev PetscInt width = xright-xleft; \ 1621eec7e3faSDmitry Karpeev *n += width*(*yhigh_loc-*ylow_loc-1); \ 1622eec7e3faSDmitry Karpeev *n += PetscMin(PetscMax(*xright_loc-xleft,0),width); \ 1623eec7e3faSDmitry Karpeev *n -= PetscMin(PetscMax(*xleft_loc-xleft,0), width); \ 1624af538404SDmitry Karpeev } \ 1625af538404SDmitry Karpeev } 1626af538404SDmitry Karpeev 1627c3518bceSDmitry Karpeev 1628eec7e3faSDmitry Karpeev 1629b1a0a8a3SJed Brown /*@ 1630f746d493SDmitry Karpeev PCGASMCreateSubdomains2D - Creates the index sets for the overlapping Schwarz 1631b1a0a8a3SJed Brown preconditioner for a two-dimensional problem on a regular grid. 1632b1a0a8a3SJed Brown 1633af538404SDmitry Karpeev Collective 1634b1a0a8a3SJed Brown 1635b1a0a8a3SJed Brown Input Parameters: 163606b43e7eSDmitry Karpeev + M, N - the global number of grid points in the x and y directions 1637af538404SDmitry Karpeev . Mdomains, Ndomains - the global number of subdomains in the x and y directions 1638b1a0a8a3SJed Brown . dof - degrees of freedom per node 1639b1a0a8a3SJed Brown - overlap - overlap in mesh lines 1640b1a0a8a3SJed Brown 1641b1a0a8a3SJed Brown Output Parameters: 1642af538404SDmitry Karpeev + Nsub - the number of local subdomains created 16436a4f0f73SDmitry Karpeev . iis - array of index sets defining inner (nonoverlapping) subdomains 16446a4f0f73SDmitry Karpeev - ois - array of index sets defining outer (overlapping, if overlap > 0) subdomains 1645b1a0a8a3SJed Brown 1646b1a0a8a3SJed Brown 1647b1a0a8a3SJed Brown Level: advanced 1648b1a0a8a3SJed Brown 1649f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, create, subdomains, 2D, regular grid 1650b1a0a8a3SJed Brown 16518f3b4b4dSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMGetSubKSP(), PCGASMSetOverlap() 1652b1a0a8a3SJed Brown @*/ 16536a4f0f73SDmitry Karpeev PetscErrorCode PCGASMCreateSubdomains2D(PC pc,PetscInt M,PetscInt N,PetscInt Mdomains,PetscInt Ndomains,PetscInt dof,PetscInt overlap,PetscInt *nsub,IS **iis,IS **ois) 1654b1a0a8a3SJed Brown { 1655b1a0a8a3SJed Brown PetscErrorCode ierr; 16562bbb417fSDmitry Karpeev PetscMPIInt size, rank; 16572bbb417fSDmitry Karpeev PetscInt i, j; 16582bbb417fSDmitry Karpeev PetscInt maxheight, maxwidth; 16592bbb417fSDmitry Karpeev PetscInt xstart, xleft, xright, xleft_loc, xright_loc; 16602bbb417fSDmitry Karpeev PetscInt ystart, ylow, yhigh, ylow_loc, yhigh_loc; 1661eec7e3faSDmitry Karpeev PetscInt x[2][2], y[2][2], n[2]; 16622bbb417fSDmitry Karpeev PetscInt first, last; 16632bbb417fSDmitry Karpeev PetscInt nidx, *idx; 16642bbb417fSDmitry Karpeev PetscInt ii,jj,s,q,d; 1665af538404SDmitry Karpeev PetscInt k,kk; 16662bbb417fSDmitry Karpeev PetscMPIInt color; 16672bbb417fSDmitry Karpeev MPI_Comm comm, subcomm; 16686a4f0f73SDmitry Karpeev IS **xis = 0, **is = ois, **is_local = iis; 1669d34fcf5fSBarry Smith 1670b1a0a8a3SJed Brown PetscFunctionBegin; 16712bbb417fSDmitry Karpeev ierr = PetscObjectGetComm((PetscObject)pc, &comm);CHKERRQ(ierr); 16722bbb417fSDmitry Karpeev ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 16732bbb417fSDmitry Karpeev ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 16742bbb417fSDmitry Karpeev ierr = MatGetOwnershipRange(pc->pmat, &first, &last);CHKERRQ(ierr); 1675d34fcf5fSBarry 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) " 16762bbb417fSDmitry Karpeev "does not respect the number of degrees of freedom per grid point %D", first, last, dof); 1677d34fcf5fSBarry Smith 1678af538404SDmitry Karpeev /* Determine the number of domains with nonzero intersections with the local ownership range. */ 16792bbb417fSDmitry Karpeev s = 0; 1680af538404SDmitry Karpeev ystart = 0; 1681af538404SDmitry Karpeev for (j=0; j<Ndomains; ++j) { 1682af538404SDmitry Karpeev maxheight = N/Ndomains + ((N % Ndomains) > j); /* Maximal height of subdomain */ 1683af538404SDmitry 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); 1684eec7e3faSDmitry Karpeev /* Vertical domain limits with an overlap. */ 1685eec7e3faSDmitry Karpeev ylow = PetscMax(ystart - overlap,0); 1686eec7e3faSDmitry Karpeev yhigh = PetscMin(ystart + maxheight + overlap,N); 1687b1a0a8a3SJed Brown xstart = 0; 1688af538404SDmitry Karpeev for (i=0; i<Mdomains; ++i) { 1689af538404SDmitry Karpeev maxwidth = M/Mdomains + ((M % Mdomains) > i); /* Maximal width of subdomain */ 1690af538404SDmitry 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); 1691eec7e3faSDmitry Karpeev /* Horizontal domain limits with an overlap. */ 1692eec7e3faSDmitry Karpeev xleft = PetscMax(xstart - overlap,0); 1693eec7e3faSDmitry Karpeev xright = PetscMin(xstart + maxwidth + overlap,M); 1694af538404SDmitry Karpeev /* 1695af538404SDmitry Karpeev Determine whether this subdomain intersects this processor's ownership range of pc->pmat. 1696af538404SDmitry Karpeev */ 1697c3518bceSDmitry Karpeev PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,(&xleft_loc),(&ylow_loc),(&xright_loc),(&yhigh_loc),(&nidx)); 16982fa5cd67SKarl Rupp if (nidx) ++s; 1699af538404SDmitry Karpeev xstart += maxwidth; 1700af538404SDmitry Karpeev } /* for (i = 0; i < Mdomains; ++i) */ 1701af538404SDmitry Karpeev ystart += maxheight; 1702af538404SDmitry Karpeev } /* for (j = 0; j < Ndomains; ++j) */ 17032fa5cd67SKarl Rupp 1704af538404SDmitry Karpeev /* Now we can allocate the necessary number of ISs. */ 1705af538404SDmitry Karpeev *nsub = s; 1706854ce69bSBarry Smith ierr = PetscMalloc1(*nsub,is);CHKERRQ(ierr); 1707854ce69bSBarry Smith ierr = PetscMalloc1(*nsub,is_local);CHKERRQ(ierr); 1708af538404SDmitry Karpeev s = 0; 1709af538404SDmitry Karpeev ystart = 0; 1710af538404SDmitry Karpeev for (j=0; j<Ndomains; ++j) { 1711af538404SDmitry Karpeev maxheight = N/Ndomains + ((N % Ndomains) > j); /* Maximal height of subdomain */ 1712af538404SDmitry 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); 1713af538404SDmitry Karpeev /* Vertical domain limits with an overlap. */ 1714eec7e3faSDmitry Karpeev y[0][0] = PetscMax(ystart - overlap,0); 1715eec7e3faSDmitry Karpeev y[0][1] = PetscMin(ystart + maxheight + overlap,N); 1716eec7e3faSDmitry Karpeev /* Vertical domain limits without an overlap. */ 1717eec7e3faSDmitry Karpeev y[1][0] = ystart; 1718eec7e3faSDmitry Karpeev y[1][1] = PetscMin(ystart + maxheight,N); 1719eec7e3faSDmitry Karpeev xstart = 0; 1720af538404SDmitry Karpeev for (i=0; i<Mdomains; ++i) { 1721af538404SDmitry Karpeev maxwidth = M/Mdomains + ((M % Mdomains) > i); /* Maximal width of subdomain */ 1722af538404SDmitry 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); 1723af538404SDmitry Karpeev /* Horizontal domain limits with an overlap. */ 1724eec7e3faSDmitry Karpeev x[0][0] = PetscMax(xstart - overlap,0); 1725eec7e3faSDmitry Karpeev x[0][1] = PetscMin(xstart + maxwidth + overlap,M); 1726eec7e3faSDmitry Karpeev /* Horizontal domain limits without an overlap. */ 1727eec7e3faSDmitry Karpeev x[1][0] = xstart; 1728eec7e3faSDmitry Karpeev x[1][1] = PetscMin(xstart+maxwidth,M); 17292bbb417fSDmitry Karpeev /* 17302bbb417fSDmitry Karpeev Determine whether this domain intersects this processor's ownership range of pc->pmat. 17312bbb417fSDmitry Karpeev Do this twice: first for the domains with overlaps, and once without. 17322bbb417fSDmitry Karpeev During the first pass create the subcommunicators, and use them on the second pass as well. 17332bbb417fSDmitry Karpeev */ 17342bbb417fSDmitry Karpeev for (q = 0; q < 2; ++q) { 1735cc96b2e5SBarry Smith PetscBool split = PETSC_FALSE; 17362bbb417fSDmitry Karpeev /* 17372bbb417fSDmitry Karpeev domain limits, (xleft, xright) and (ylow, yheigh) are adjusted 17382bbb417fSDmitry Karpeev according to whether the domain with an overlap or without is considered. 17392bbb417fSDmitry Karpeev */ 17402bbb417fSDmitry Karpeev xleft = x[q][0]; xright = x[q][1]; 17412bbb417fSDmitry Karpeev ylow = y[q][0]; yhigh = y[q][1]; 1742c3518bceSDmitry Karpeev PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,(&xleft_loc),(&ylow_loc),(&xright_loc),(&yhigh_loc),(&nidx)); 1743af538404SDmitry Karpeev nidx *= dof; 1744eec7e3faSDmitry Karpeev n[q] = nidx; 17452bbb417fSDmitry Karpeev /* 1746eec7e3faSDmitry Karpeev Based on the counted number of indices in the local domain *with an overlap*, 1747af538404SDmitry Karpeev construct a subcommunicator of all the processors supporting this domain. 1748eec7e3faSDmitry Karpeev Observe that a domain with an overlap might have nontrivial local support, 1749eec7e3faSDmitry Karpeev while the domain without an overlap might not. Hence, the decision to participate 1750eec7e3faSDmitry Karpeev in the subcommunicator must be based on the domain with an overlap. 17512bbb417fSDmitry Karpeev */ 17522bbb417fSDmitry Karpeev if (q == 0) { 17532fa5cd67SKarl Rupp if (nidx) color = 1; 17542fa5cd67SKarl Rupp else color = MPI_UNDEFINED; 17552bbb417fSDmitry Karpeev ierr = MPI_Comm_split(comm, color, rank, &subcomm);CHKERRQ(ierr); 1756cc96b2e5SBarry Smith split = PETSC_TRUE; 17572bbb417fSDmitry Karpeev } 1758af538404SDmitry Karpeev /* 1759eec7e3faSDmitry Karpeev Proceed only if the number of local indices *with an overlap* is nonzero. 1760af538404SDmitry Karpeev */ 1761eec7e3faSDmitry Karpeev if (n[0]) { 17622fa5cd67SKarl Rupp if (q == 0) xis = is; 1763af538404SDmitry Karpeev if (q == 1) { 1764af538404SDmitry Karpeev /* 1765eec7e3faSDmitry Karpeev The IS for the no-overlap subdomain shares a communicator with the overlapping domain. 1766af538404SDmitry Karpeev Moreover, if the overlap is zero, the two ISs are identical. 1767af538404SDmitry Karpeev */ 1768b1a0a8a3SJed Brown if (overlap == 0) { 1769eec7e3faSDmitry Karpeev (*is_local)[s] = (*is)[s]; 17702bbb417fSDmitry Karpeev ierr = PetscObjectReference((PetscObject)(*is)[s]);CHKERRQ(ierr); 17712bbb417fSDmitry Karpeev continue; 1772d34fcf5fSBarry Smith } else { 17736a4f0f73SDmitry Karpeev xis = is_local; 1774eec7e3faSDmitry Karpeev subcomm = ((PetscObject)(*is)[s])->comm; 17752bbb417fSDmitry Karpeev } 1776af538404SDmitry Karpeev } /* if (q == 1) */ 17770298fd71SBarry Smith idx = NULL; 1778785e854fSJed Brown ierr = PetscMalloc1(nidx,&idx);CHKERRQ(ierr); 1779eec7e3faSDmitry Karpeev if (nidx) { 17802bbb417fSDmitry Karpeev k = 0; 17812bbb417fSDmitry Karpeev for (jj=ylow_loc; jj<yhigh_loc; ++jj) { 1782af538404SDmitry Karpeev PetscInt x0 = (jj==ylow_loc) ? xleft_loc : xleft; 1783af538404SDmitry Karpeev PetscInt x1 = (jj==yhigh_loc-1) ? xright_loc : xright; 1784af538404SDmitry Karpeev kk = dof*(M*jj + x0); 17852bbb417fSDmitry Karpeev for (ii=x0; ii<x1; ++ii) { 17862bbb417fSDmitry Karpeev for (d = 0; d < dof; ++d) { 17872bbb417fSDmitry Karpeev idx[k++] = kk++; 1788b1a0a8a3SJed Brown } 1789b1a0a8a3SJed Brown } 1790b1a0a8a3SJed Brown } 1791eec7e3faSDmitry Karpeev } 17926a4f0f73SDmitry Karpeev ierr = ISCreateGeneral(subcomm,nidx,idx,PETSC_OWN_POINTER,(*xis)+s);CHKERRQ(ierr); 1793cc96b2e5SBarry Smith if (split) { 1794cc96b2e5SBarry Smith ierr = MPI_Comm_free(&subcomm);CHKERRQ(ierr); 1795cc96b2e5SBarry Smith } 1796eec7e3faSDmitry Karpeev }/* if (n[0]) */ 17972bbb417fSDmitry Karpeev }/* for (q = 0; q < 2; ++q) */ 17982fa5cd67SKarl Rupp if (n[0]) ++s; 1799af538404SDmitry Karpeev xstart += maxwidth; 1800af538404SDmitry Karpeev } /* for (i = 0; i < Mdomains; ++i) */ 18012bbb417fSDmitry Karpeev ystart += maxheight; 1802af538404SDmitry Karpeev } /* for (j = 0; j < Ndomains; ++j) */ 1803b1a0a8a3SJed Brown PetscFunctionReturn(0); 1804b1a0a8a3SJed Brown } 1805b1a0a8a3SJed Brown 1806b1a0a8a3SJed Brown /*@C 180706b43e7eSDmitry Karpeev PCGASMGetSubdomains - Gets the subdomains supported on this processor 180806b43e7eSDmitry Karpeev for the additive Schwarz preconditioner. 1809b1a0a8a3SJed Brown 1810b1a0a8a3SJed Brown Not Collective 1811b1a0a8a3SJed Brown 1812b1a0a8a3SJed Brown Input Parameter: 1813b1a0a8a3SJed Brown . pc - the preconditioner context 1814b1a0a8a3SJed Brown 1815b1a0a8a3SJed Brown Output Parameters: 1816b1a0a8a3SJed Brown + n - the number of subdomains for this processor (default value = 1) 18170298fd71SBarry Smith . iis - the index sets that define the inner subdomains (without overlap) supported on this processor (can be NULL) 18180298fd71SBarry Smith - ois - the index sets that define the outer subdomains (with overlap) supported on this processor (can be NULL) 1819b1a0a8a3SJed Brown 1820b1a0a8a3SJed Brown 1821b1a0a8a3SJed Brown Notes: 18226a4f0f73SDmitry Karpeev The user is responsible for destroying the ISs and freeing the returned arrays. 1823b1a0a8a3SJed Brown The IS numbering is in the parallel, global numbering of the vector. 1824b1a0a8a3SJed Brown 1825b1a0a8a3SJed Brown Level: advanced 1826b1a0a8a3SJed Brown 182706b43e7eSDmitry Karpeev .keywords: PC, GASM, get, subdomains, additive Schwarz 1828b1a0a8a3SJed Brown 18298f3b4b4dSDmitry Karpeev .seealso: PCGASMSetOverlap(), PCGASMGetSubKSP(), PCGASMCreateSubdomains2D(), 18308f3b4b4dSDmitry Karpeev PCGASMSetSubdomains(), PCGASMGetSubmatrices() 1831b1a0a8a3SJed Brown @*/ 18326a4f0f73SDmitry Karpeev PetscErrorCode PCGASMGetSubdomains(PC pc,PetscInt *n,IS *iis[],IS *ois[]) 1833b1a0a8a3SJed Brown { 1834f746d493SDmitry Karpeev PC_GASM *osm; 1835b1a0a8a3SJed Brown PetscErrorCode ierr; 1836b1a0a8a3SJed Brown PetscBool match; 18376a4f0f73SDmitry Karpeev PetscInt i; 18385fd66863SKarl Rupp 1839b1a0a8a3SJed Brown PetscFunctionBegin; 1840b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1841251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match);CHKERRQ(ierr); 1842ce94432eSBarry Smith if (!match) SETERRQ2(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Incorrect object type: expected %s, got %s instead", PCGASM, ((PetscObject)pc)->type_name); 1843f746d493SDmitry Karpeev osm = (PC_GASM*)pc->data; 1844ab3e923fSDmitry Karpeev if (n) *n = osm->n; 18456a4f0f73SDmitry Karpeev if (iis) { 1846785e854fSJed Brown ierr = PetscMalloc1(osm->n, iis);CHKERRQ(ierr); 18476a4f0f73SDmitry Karpeev } 18486a4f0f73SDmitry Karpeev if (ois) { 1849785e854fSJed Brown ierr = PetscMalloc1(osm->n, ois);CHKERRQ(ierr); 18506a4f0f73SDmitry Karpeev } 18516a4f0f73SDmitry Karpeev if (iis || ois) { 18526a4f0f73SDmitry Karpeev for (i = 0; i < osm->n; ++i) { 18536a4f0f73SDmitry Karpeev if (iis) (*iis)[i] = osm->iis[i]; 18546a4f0f73SDmitry Karpeev if (ois) (*ois)[i] = osm->ois[i]; 18556a4f0f73SDmitry Karpeev } 1856b1a0a8a3SJed Brown } 1857b1a0a8a3SJed Brown PetscFunctionReturn(0); 1858b1a0a8a3SJed Brown } 1859b1a0a8a3SJed Brown 1860b1a0a8a3SJed Brown /*@C 186106b43e7eSDmitry Karpeev PCGASMGetSubmatrices - Gets the local submatrices (for this processor 1862b1a0a8a3SJed Brown only) for the additive Schwarz preconditioner. 1863b1a0a8a3SJed Brown 1864b1a0a8a3SJed Brown Not Collective 1865b1a0a8a3SJed Brown 1866b1a0a8a3SJed Brown Input Parameter: 1867b1a0a8a3SJed Brown . pc - the preconditioner context 1868b1a0a8a3SJed Brown 1869b1a0a8a3SJed Brown Output Parameters: 1870b1a0a8a3SJed Brown + n - the number of matrices for this processor (default value = 1) 1871b1a0a8a3SJed Brown - mat - the matrices 1872b1a0a8a3SJed Brown 187395452b02SPatrick Sanan Notes: 187495452b02SPatrick Sanan matrices returned by this routine have the same communicators as the index sets (IS) 18758f3b4b4dSDmitry Karpeev used to define subdomains in PCGASMSetSubdomains() 1876b1a0a8a3SJed Brown Level: advanced 1877b1a0a8a3SJed Brown 1878f746d493SDmitry Karpeev .keywords: PC, GASM, set, local, subdomains, additive Schwarz, block Jacobi 1879b1a0a8a3SJed Brown 18808f3b4b4dSDmitry Karpeev .seealso: PCGASMSetOverlap(), PCGASMGetSubKSP(), 188106b43e7eSDmitry Karpeev PCGASMCreateSubdomains2D(), PCGASMSetSubdomains(), PCGASMGetSubdomains() 1882b1a0a8a3SJed Brown @*/ 188306b43e7eSDmitry Karpeev PetscErrorCode PCGASMGetSubmatrices(PC pc,PetscInt *n,Mat *mat[]) 1884b1a0a8a3SJed Brown { 1885f746d493SDmitry Karpeev PC_GASM *osm; 1886b1a0a8a3SJed Brown PetscErrorCode ierr; 1887b1a0a8a3SJed Brown PetscBool match; 1888b1a0a8a3SJed Brown 1889b1a0a8a3SJed Brown PetscFunctionBegin; 1890b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1891b1a0a8a3SJed Brown PetscValidIntPointer(n,2); 1892b1a0a8a3SJed Brown if (mat) PetscValidPointer(mat,3); 1893ce94432eSBarry Smith if (!pc->setupcalled) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call after KSPSetUP() or PCSetUp()."); 1894251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match);CHKERRQ(ierr); 1895ce94432eSBarry Smith if (!match) SETERRQ2(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Expected %s, got %s instead", PCGASM, ((PetscObject)pc)->type_name); 1896f746d493SDmitry Karpeev osm = (PC_GASM*)pc->data; 1897ab3e923fSDmitry Karpeev if (n) *n = osm->n; 1898b1a0a8a3SJed Brown if (mat) *mat = osm->pmat; 1899b1a0a8a3SJed Brown PetscFunctionReturn(0); 1900b1a0a8a3SJed Brown } 1901d709ea83SDmitry Karpeev 1902d709ea83SDmitry Karpeev /*@ 19038f3b4b4dSDmitry Karpeev PCGASMSetUseDMSubdomains - Indicates whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible. 1904d709ea83SDmitry Karpeev Logically Collective 1905d709ea83SDmitry Karpeev 1906d709ea83SDmitry Karpeev Input Parameter: 1907d709ea83SDmitry Karpeev + pc - the preconditioner 1908d709ea83SDmitry Karpeev - flg - boolean indicating whether to use subdomains defined by the DM 1909d709ea83SDmitry Karpeev 1910d709ea83SDmitry Karpeev Options Database Key: 19118f3b4b4dSDmitry Karpeev . -pc_gasm_dm_subdomains -pc_gasm_overlap -pc_gasm_total_subdomains 1912d709ea83SDmitry Karpeev 1913d709ea83SDmitry Karpeev Level: intermediate 1914d709ea83SDmitry Karpeev 1915d709ea83SDmitry Karpeev Notes: 19168f3b4b4dSDmitry Karpeev PCGASMSetSubdomains(), PCGASMSetTotalSubdomains() or PCGASMSetOverlap() take precedence over PCGASMSetUseDMSubdomains(), 19178f3b4b4dSDmitry Karpeev so setting PCGASMSetSubdomains() with nontrivial subdomain ISs or any of PCGASMSetTotalSubdomains() and PCGASMSetOverlap() 19188f3b4b4dSDmitry Karpeev automatically turns the latter off. 1919d709ea83SDmitry Karpeev 1920d709ea83SDmitry Karpeev .keywords: PC, ASM, DM, set, subdomains, additive Schwarz 1921d709ea83SDmitry Karpeev 19228f3b4b4dSDmitry Karpeev .seealso: PCGASMGetUseDMSubdomains(), PCGASMSetSubdomains(), PCGASMSetOverlap() 1923d709ea83SDmitry Karpeev PCGASMCreateSubdomains2D() 1924d709ea83SDmitry Karpeev @*/ 19258f3b4b4dSDmitry Karpeev PetscErrorCode PCGASMSetUseDMSubdomains(PC pc,PetscBool flg) 1926d709ea83SDmitry Karpeev { 1927d709ea83SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 1928d709ea83SDmitry Karpeev PetscErrorCode ierr; 1929d709ea83SDmitry Karpeev PetscBool match; 1930d709ea83SDmitry Karpeev 1931d709ea83SDmitry Karpeev PetscFunctionBegin; 1932d709ea83SDmitry Karpeev PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1933d709ea83SDmitry Karpeev PetscValidLogicalCollectiveBool(pc,flg,2); 1934d709ea83SDmitry Karpeev if (pc->setupcalled) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Not for a setup PC."); 1935d709ea83SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match);CHKERRQ(ierr); 1936d709ea83SDmitry Karpeev if (match) { 19378f3b4b4dSDmitry Karpeev if (!osm->user_subdomains && osm->N == PETSC_DETERMINE && osm->overlap < 0) { 1938d709ea83SDmitry Karpeev osm->dm_subdomains = flg; 1939d709ea83SDmitry Karpeev } 19408f3b4b4dSDmitry Karpeev } 1941d709ea83SDmitry Karpeev PetscFunctionReturn(0); 1942d709ea83SDmitry Karpeev } 1943d709ea83SDmitry Karpeev 1944d709ea83SDmitry Karpeev /*@ 19458f3b4b4dSDmitry Karpeev PCGASMGetUseDMSubdomains - Returns flag indicating whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible. 1946d709ea83SDmitry Karpeev Not Collective 1947d709ea83SDmitry Karpeev 1948d709ea83SDmitry Karpeev Input Parameter: 1949d709ea83SDmitry Karpeev . pc - the preconditioner 1950d709ea83SDmitry Karpeev 1951d709ea83SDmitry Karpeev Output Parameter: 1952d709ea83SDmitry Karpeev . flg - boolean indicating whether to use subdomains defined by the DM 1953d709ea83SDmitry Karpeev 1954d709ea83SDmitry Karpeev Level: intermediate 1955d709ea83SDmitry Karpeev 1956d709ea83SDmitry Karpeev .keywords: PC, ASM, DM, set, subdomains, additive Schwarz 1957d709ea83SDmitry Karpeev 19588f3b4b4dSDmitry Karpeev .seealso: PCGASMSetUseDMSubdomains(), PCGASMSetOverlap() 1959d709ea83SDmitry Karpeev PCGASMCreateSubdomains2D() 1960d709ea83SDmitry Karpeev @*/ 19618f3b4b4dSDmitry Karpeev PetscErrorCode PCGASMGetUseDMSubdomains(PC pc,PetscBool* flg) 1962d709ea83SDmitry Karpeev { 1963d709ea83SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 1964d709ea83SDmitry Karpeev PetscErrorCode ierr; 1965d709ea83SDmitry Karpeev PetscBool match; 1966d709ea83SDmitry Karpeev 1967d709ea83SDmitry Karpeev PetscFunctionBegin; 1968d709ea83SDmitry Karpeev PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1969d709ea83SDmitry Karpeev PetscValidPointer(flg,2); 1970d709ea83SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match);CHKERRQ(ierr); 1971d709ea83SDmitry Karpeev if (match) { 1972d709ea83SDmitry Karpeev if (flg) *flg = osm->dm_subdomains; 1973d709ea83SDmitry Karpeev } 1974d709ea83SDmitry Karpeev PetscFunctionReturn(0); 1975d709ea83SDmitry Karpeev } 1976