1b1a0a8a3SJed Brown /* 2f746d493SDmitry Karpeev This file defines an "generalized" additive Schwarz preconditioner for any Mat implementation. 36a4f0f73SDmitry Karpeev In this version each processor may intersect multiple subdomains and any subdomain may 46a4f0f73SDmitry Karpeev intersect multiple processors. Intersections of subdomains with processors are called *local 56a4f0f73SDmitry Karpeev subdomains*. 6b1a0a8a3SJed Brown 78f3b4b4dSDmitry Karpeev N - total number of distinct global subdomains (set explicitly in PCGASMSetTotalSubdomains() or implicitly PCGASMSetSubdomains() and then calculated in PCSetUp_GASM()) 86a4f0f73SDmitry Karpeev n - actual number of local subdomains on this processor (set in PCGASMSetSubdomains() or calculated in PCGASMSetTotalSubdomains()) 98f3b4b4dSDmitry Karpeev nmax - maximum number of local subdomains per processor (calculated in PCSetUp_GASM()) 10b1a0a8a3SJed Brown */ 11af0996ceSBarry Smith #include <petsc/private/pcimpl.h> /*I "petscpc.h" I*/ 121e25c274SJed Brown #include <petscdm.h> 13b1a0a8a3SJed Brown 14b1a0a8a3SJed Brown typedef struct { 15f746d493SDmitry Karpeev PetscInt N,n,nmax; 16b1a0a8a3SJed Brown PetscInt overlap; /* overlap requested by user */ 178f3b4b4dSDmitry Karpeev PCGASMType type; /* use reduced interpolation, restriction or both */ 188f3b4b4dSDmitry Karpeev PetscBool type_set; /* if user set this value (so won't change it for symmetric problems) */ 198f3b4b4dSDmitry Karpeev PetscBool same_subdomain_solvers; /* flag indicating whether all local solvers are same */ 208f3b4b4dSDmitry Karpeev PetscBool sort_indices; /* flag to sort subdomain indices */ 218f3b4b4dSDmitry Karpeev PetscBool user_subdomains; /* whether the user set explicit subdomain index sets -- keep them on PCReset() */ 228f3b4b4dSDmitry Karpeev PetscBool dm_subdomains; /* whether DM is allowed to define subdomains */ 238f3b4b4dSDmitry Karpeev IS *ois; /* index sets that define the outer (conceptually, overlapping) subdomains */ 248f3b4b4dSDmitry Karpeev IS *iis; /* index sets that define the inner (conceptually, nonoverlapping) subdomains */ 258f3b4b4dSDmitry Karpeev KSP *ksp; /* linear solvers for each subdomain */ 268f3b4b4dSDmitry Karpeev Mat *pmat; /* subdomain block matrices */ 27f746d493SDmitry Karpeev Vec gx,gy; /* Merged work vectors */ 28f746d493SDmitry Karpeev Vec *x,*y; /* Split work vectors; storage aliases pieces of storage of the above merged vectors. */ 296a4f0f73SDmitry Karpeev VecScatter gorestriction; /* merged restriction to disjoint union of outer subdomains */ 306a4f0f73SDmitry Karpeev VecScatter girestriction; /* merged restriction to disjoint union of inner subdomains */ 31f746d493SDmitry Karpeev } PC_GASM; 32b1a0a8a3SJed Brown 33b1a0a8a3SJed Brown #undef __FUNCT__ 348f3b4b4dSDmitry Karpeev #define __FUNCT__ "PCGASMComputeGlobalSubdomainNumbering_Private" 358f3b4b4dSDmitry Karpeev static PetscErrorCode PCGASMComputeGlobalSubdomainNumbering_Private(PC pc,PetscInt **numbering,PetscInt **permutation) 368f3b4b4dSDmitry Karpeev { 378f3b4b4dSDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 388f3b4b4dSDmitry Karpeev PetscInt i; 398f3b4b4dSDmitry Karpeev PetscErrorCode ierr; 408f3b4b4dSDmitry Karpeev 418f3b4b4dSDmitry Karpeev PetscFunctionBegin; 428f3b4b4dSDmitry Karpeev /* Determine the number of globally-distinct subdomains and compute a global numbering for them. */ 438f3b4b4dSDmitry Karpeev ierr = PetscMalloc2(osm->n,numbering,osm->n,permutation);CHKERRQ(ierr); 448f3b4b4dSDmitry Karpeev ierr = PetscObjectsListGetGlobalNumbering(PetscObjectComm((PetscObject)pc),osm->n,(PetscObject*)osm->iis,NULL,*numbering);CHKERRQ(ierr); 458f3b4b4dSDmitry Karpeev for (i = 0; i < osm->n; ++i) (*permutation)[i] = i; 468f3b4b4dSDmitry Karpeev ierr = PetscSortIntWithPermutation(osm->n,*numbering,*permutation);CHKERRQ(ierr); 478f3b4b4dSDmitry Karpeev PetscFunctionReturn(0); 488f3b4b4dSDmitry Karpeev } 498f3b4b4dSDmitry Karpeev 508f3b4b4dSDmitry Karpeev #undef __FUNCT__ 5106b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMSubdomainView_Private" 5206b43e7eSDmitry Karpeev static PetscErrorCode PCGASMSubdomainView_Private(PC pc, PetscInt i, PetscViewer viewer) 53af538404SDmitry Karpeev { 54af538404SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 5506b43e7eSDmitry Karpeev PetscInt j,nidx; 56af538404SDmitry Karpeev const PetscInt *idx; 5706b43e7eSDmitry Karpeev PetscViewer sviewer; 58af538404SDmitry Karpeev char *cidx; 59af538404SDmitry Karpeev PetscErrorCode ierr; 60af538404SDmitry Karpeev 61af538404SDmitry Karpeev PetscFunctionBegin; 62ce94432eSBarry Smith if (i < 0 || i > osm->n) SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONG, "Invalid subdomain %D: must nonnegative and less than %D", i, osm->n); 634bde246dSDmitry Karpeev /* Inner subdomains. */ 644bde246dSDmitry Karpeev ierr = ISGetLocalSize(osm->iis[i], &nidx);CHKERRQ(ierr); 654bde246dSDmitry Karpeev /* 664bde246dSDmitry Karpeev No more than 15 characters per index plus a space. 674bde246dSDmitry Karpeev PetscViewerStringSPrintf requires a string of size at least 2, so use (nidx+1) instead of nidx, 684bde246dSDmitry Karpeev in case nidx == 0. That will take care of the space for the trailing '\0' as well. 694bde246dSDmitry Karpeev For nidx == 0, the whole string 16 '\0'. 704bde246dSDmitry Karpeev */ 7189d949e2SBarry Smith ierr = PetscMalloc1(16*(nidx+1)+1, &cidx);CHKERRQ(ierr); 724bde246dSDmitry Karpeev ierr = PetscViewerStringOpen(PETSC_COMM_SELF, cidx, 16*(nidx+1)+1, &sviewer);CHKERRQ(ierr); 734bde246dSDmitry Karpeev ierr = ISGetIndices(osm->iis[i], &idx);CHKERRQ(ierr); 744bde246dSDmitry Karpeev for (j = 0; j < nidx; ++j) { 754bde246dSDmitry Karpeev ierr = PetscViewerStringSPrintf(sviewer, "%D ", idx[j]);CHKERRQ(ierr); 764bde246dSDmitry Karpeev } 774bde246dSDmitry Karpeev ierr = ISRestoreIndices(osm->iis[i],&idx);CHKERRQ(ierr); 784bde246dSDmitry Karpeev ierr = PetscViewerDestroy(&sviewer);CHKERRQ(ierr); 794bde246dSDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer, "Inner subdomain:\n");CHKERRQ(ierr); 804bde246dSDmitry Karpeev ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 81*c5e4d11fSDmitry Karpeev ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 824bde246dSDmitry Karpeev ierr = PetscViewerASCIISynchronizedPrintf(viewer, "%s", cidx);CHKERRQ(ierr); 834bde246dSDmitry Karpeev ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 84*c5e4d11fSDmitry Karpeev ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 854bde246dSDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer, "\n");CHKERRQ(ierr); 864bde246dSDmitry Karpeev ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 874bde246dSDmitry Karpeev ierr = PetscFree(cidx);CHKERRQ(ierr); 884bde246dSDmitry Karpeev /* Outer subdomains. */ 896a4f0f73SDmitry Karpeev ierr = ISGetLocalSize(osm->ois[i], &nidx);CHKERRQ(ierr); 90eec7e3faSDmitry Karpeev /* 91eec7e3faSDmitry Karpeev No more than 15 characters per index plus a space. 92eec7e3faSDmitry Karpeev PetscViewerStringSPrintf requires a string of size at least 2, so use (nidx+1) instead of nidx, 93eec7e3faSDmitry Karpeev in case nidx == 0. That will take care of the space for the trailing '\0' as well. 94eec7e3faSDmitry Karpeev For nidx == 0, the whole string 16 '\0'. 95eec7e3faSDmitry Karpeev */ 9689d949e2SBarry Smith ierr = PetscMalloc1(16*(nidx+1)+1, &cidx);CHKERRQ(ierr); 9706b43e7eSDmitry Karpeev ierr = PetscViewerStringOpen(PETSC_COMM_SELF, cidx, 16*(nidx+1)+1, &sviewer);CHKERRQ(ierr); 986a4f0f73SDmitry Karpeev ierr = ISGetIndices(osm->ois[i], &idx);CHKERRQ(ierr); 99af538404SDmitry Karpeev for (j = 0; j < nidx; ++j) { 100af538404SDmitry Karpeev ierr = PetscViewerStringSPrintf(sviewer,"%D ", idx[j]);CHKERRQ(ierr); 101af538404SDmitry Karpeev } 1026bf464f9SBarry Smith ierr = PetscViewerDestroy(&sviewer);CHKERRQ(ierr); 1036a4f0f73SDmitry Karpeev ierr = ISRestoreIndices(osm->ois[i],&idx);CHKERRQ(ierr); 1046a4f0f73SDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer, "Outer subdomain:\n");CHKERRQ(ierr); 10506b43e7eSDmitry Karpeev ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 106*c5e4d11fSDmitry Karpeev ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 107af538404SDmitry Karpeev ierr = PetscViewerASCIISynchronizedPrintf(viewer, "%s", cidx);CHKERRQ(ierr); 108af538404SDmitry Karpeev ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 109*c5e4d11fSDmitry Karpeev ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 110af538404SDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer, "\n");CHKERRQ(ierr); 11106b43e7eSDmitry Karpeev ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 112af538404SDmitry Karpeev ierr = PetscFree(cidx);CHKERRQ(ierr); 11306b43e7eSDmitry Karpeev PetscFunctionReturn(0); 11406b43e7eSDmitry Karpeev } 11506b43e7eSDmitry Karpeev 11606b43e7eSDmitry Karpeev #undef __FUNCT__ 11706b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMPrintSubdomains" 11806b43e7eSDmitry Karpeev static PetscErrorCode PCGASMPrintSubdomains(PC pc) 11906b43e7eSDmitry Karpeev { 12006b43e7eSDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 12106b43e7eSDmitry Karpeev const char *prefix; 12206b43e7eSDmitry Karpeev char fname[PETSC_MAX_PATH_LEN+1]; 1238f3b4b4dSDmitry Karpeev PetscInt l, d, count; 1248f3b4b4dSDmitry Karpeev PetscBool doprint,found; 1250298fd71SBarry Smith PetscViewer viewer, sviewer = NULL; 1268f3b4b4dSDmitry Karpeev PetscInt *numbering,*permutation;/* global numbering of locally-supported subdomains and the permutation from the local ordering */ 12706b43e7eSDmitry Karpeev PetscErrorCode ierr; 12806b43e7eSDmitry Karpeev 12906b43e7eSDmitry Karpeev PetscFunctionBegin; 13006b43e7eSDmitry Karpeev ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 1318f3b4b4dSDmitry Karpeev doprint = PETSC_FALSE; 1328f3b4b4dSDmitry Karpeev ierr = PetscOptionsGetBool(prefix,"-pc_gasm_print_subdomains",&doprint,NULL);CHKERRQ(ierr); 1338f3b4b4dSDmitry Karpeev if (!doprint) PetscFunctionReturn(0); 13406b43e7eSDmitry Karpeev ierr = PetscOptionsGetString(prefix,"-pc_gasm_print_subdomains",fname,PETSC_MAX_PATH_LEN,&found);CHKERRQ(ierr); 13506b43e7eSDmitry Karpeev if (!found) { ierr = PetscStrcpy(fname,"stdout");CHKERRQ(ierr); }; 136ce94432eSBarry Smith ierr = PetscViewerASCIIOpen(PetscObjectComm((PetscObject)pc),fname,&viewer);CHKERRQ(ierr); 1374bde246dSDmitry Karpeev /* 1384bde246dSDmitry Karpeev Make sure the viewer has a name. Otherwise this may cause a deadlock or other weird errors when creating a subcomm viewer: 1394bde246dSDmitry Karpeev the subcomm viewer will attempt to inherit the viewer's name, which, if not set, will be constructed collectively on the comm. 1404bde246dSDmitry Karpeev */ 1414bde246dSDmitry Karpeev ierr = PetscObjectName((PetscObject)viewer);CHKERRQ(ierr); 1424bde246dSDmitry Karpeev l = 0; 1438f3b4b4dSDmitry Karpeev ierr = PCGASMComputeGlobalSubdomainNumbering_Private(pc,&numbering,&permutation);CHKERRQ(ierr); 1448f3b4b4dSDmitry Karpeev for (count = 0; count < osm->N; ++count) { 1454bde246dSDmitry Karpeev /* Now let subdomains go one at a time in the global numbering order and print their subdomain/solver info. */ 1464bde246dSDmitry Karpeev if (l<osm->n) { 1474bde246dSDmitry Karpeev d = permutation[l]; /* d is the local number of the l-th smallest (in the global ordering) among the locally supported subdomains */ 1484bde246dSDmitry Karpeev if (numbering[d] == count) { 149*c5e4d11fSDmitry Karpeev ierr = PetscViewerGetSubViewer(viewer,((PetscObject)osm->ois[d])->comm, &sviewer);CHKERRQ(ierr); 1504bde246dSDmitry Karpeev ierr = PCGASMSubdomainView_Private(pc,d,sviewer);CHKERRQ(ierr); 151*c5e4d11fSDmitry Karpeev ierr = PetscViewerRestoreSubViewer(viewer,((PetscObject)osm->ois[d])->comm, &sviewer);CHKERRQ(ierr); 1524bde246dSDmitry Karpeev ++l; 153af538404SDmitry Karpeev } 1544bde246dSDmitry Karpeev } 155ce94432eSBarry Smith ierr = MPI_Barrier(PetscObjectComm((PetscObject)pc));CHKERRQ(ierr); 1564bde246dSDmitry Karpeev } 1578f3b4b4dSDmitry Karpeev ierr = PetscFree2(numbering,permutation);CHKERRQ(ierr); 1584bde246dSDmitry Karpeev ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 159af538404SDmitry Karpeev PetscFunctionReturn(0); 160af538404SDmitry Karpeev } 161af538404SDmitry Karpeev 162af538404SDmitry Karpeev 163af538404SDmitry Karpeev #undef __FUNCT__ 164f746d493SDmitry Karpeev #define __FUNCT__ "PCView_GASM" 165f746d493SDmitry Karpeev static PetscErrorCode PCView_GASM(PC pc,PetscViewer viewer) 166b1a0a8a3SJed Brown { 167f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 168af538404SDmitry Karpeev const char *prefix; 169b1a0a8a3SJed Brown PetscErrorCode ierr; 170af538404SDmitry Karpeev PetscMPIInt rank, size; 1718f3b4b4dSDmitry Karpeev PetscInt bsz; 17206b43e7eSDmitry Karpeev PetscBool iascii,view_subdomains=PETSC_FALSE; 173b1a0a8a3SJed Brown PetscViewer sviewer; 1748f3b4b4dSDmitry Karpeev PetscInt count, l; 1756a4f0f73SDmitry Karpeev char overlap[256] = "user-defined overlap"; 1766a4f0f73SDmitry Karpeev char gsubdomains[256] = "unknown total number of subdomains"; 17706b43e7eSDmitry Karpeev char msubdomains[256] = "unknown max number of local subdomains"; 1788f3b4b4dSDmitry Karpeev PetscInt *numbering,*permutation;/* global numbering of locally-supported subdomains and the permutation from the local ordering */ 17911aeaf0aSBarry Smith 180b1a0a8a3SJed Brown PetscFunctionBegin; 181ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size);CHKERRQ(ierr); 182ce94432eSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank);CHKERRQ(ierr); 18306b43e7eSDmitry Karpeev 18406b43e7eSDmitry Karpeev if (osm->overlap >= 0) { 18506b43e7eSDmitry Karpeev ierr = PetscSNPrintf(overlap,sizeof(overlap),"requested amount of overlap = %D",osm->overlap);CHKERRQ(ierr); 18606b43e7eSDmitry Karpeev } 1878f3b4b4dSDmitry Karpeev if (osm->N != PETSC_DETERMINE) { 1888f3b4b4dSDmitry Karpeev ierr = PetscSNPrintf(gsubdomains, sizeof(gsubdomains), "total number of subdomains = %D",osm->N);CHKERRQ(ierr); 18906b43e7eSDmitry Karpeev } 1908f3b4b4dSDmitry Karpeev if (osm->nmax != PETSC_DETERMINE) { 19106b43e7eSDmitry Karpeev ierr = PetscSNPrintf(msubdomains,sizeof(msubdomains),"max number of local subdomains = %D",osm->nmax);CHKERRQ(ierr); 19206b43e7eSDmitry Karpeev } 19306b43e7eSDmitry Karpeev 194af538404SDmitry Karpeev ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 1950298fd71SBarry Smith ierr = PetscOptionsGetBool(prefix,"-pc_gasm_view_subdomains",&view_subdomains,NULL);CHKERRQ(ierr); 19606b43e7eSDmitry Karpeev 19706b43e7eSDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 198b1a0a8a3SJed Brown if (iascii) { 19906b43e7eSDmitry Karpeev /* 20006b43e7eSDmitry Karpeev Make sure the viewer has a name. Otherwise this may cause a deadlock when creating a subcomm viewer: 20106b43e7eSDmitry Karpeev the subcomm viewer will attempt to inherit the viewer's name, which, if not set, will be constructed 20206b43e7eSDmitry Karpeev collectively on the comm. 20306b43e7eSDmitry Karpeev */ 20406b43e7eSDmitry Karpeev ierr = PetscObjectName((PetscObject)viewer);CHKERRQ(ierr); 20506b43e7eSDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer,"Generalized additive Schwarz:\n");CHKERRQ(ierr); 20606b43e7eSDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer,"Restriction/interpolation type: %s\n",PCGASMTypes[osm->type]);CHKERRQ(ierr); 20706b43e7eSDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer,"%s\n",overlap);CHKERRQ(ierr); 20806b43e7eSDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer,"%s\n",gsubdomains);CHKERRQ(ierr); 20906b43e7eSDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer,"%s\n",msubdomains);CHKERRQ(ierr); 210*c5e4d11fSDmitry Karpeev ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 2118f3b4b4dSDmitry Karpeev ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d|%d] number of locally-supported subdomains = %D\n",rank,size,osm->n);CHKERRQ(ierr); 212af538404SDmitry Karpeev ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 213*c5e4d11fSDmitry Karpeev ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 2146a4f0f73SDmitry Karpeev /* Cannot take advantage of osm->same_subdomain_solvers without a global numbering of subdomains. */ 21506b43e7eSDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer,"Subdomain solver info is as follows:\n");CHKERRQ(ierr); 216b1a0a8a3SJed Brown ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 217b1a0a8a3SJed Brown ierr = PetscViewerASCIIPrintf(viewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr); 21806b43e7eSDmitry Karpeev /* Make sure that everybody waits for the banner to be printed. */ 219ce94432eSBarry Smith ierr = MPI_Barrier(PetscObjectComm((PetscObject)viewer));CHKERRQ(ierr); 2204bde246dSDmitry Karpeev /* Now let subdomains go one at a time in the global numbering order and print their subdomain/solver info. */ 2218f3b4b4dSDmitry Karpeev ierr = PCGASMComputeGlobalSubdomainNumbering_Private(pc,&numbering,&permutation);CHKERRQ(ierr); 22206b43e7eSDmitry Karpeev l = 0; 2238f3b4b4dSDmitry Karpeev for (count = 0; count < osm->N; ++count) { 22406b43e7eSDmitry Karpeev PetscMPIInt srank, ssize; 22506b43e7eSDmitry Karpeev if (l<osm->n) { 22606b43e7eSDmitry Karpeev PetscInt d = permutation[l]; /* d is the local number of the l-th smallest (in the global ordering) among the locally supported subdomains */ 22706b43e7eSDmitry Karpeev if (numbering[d] == count) { 2286a4f0f73SDmitry Karpeev ierr = MPI_Comm_size(((PetscObject)osm->ois[d])->comm, &ssize);CHKERRQ(ierr); 2296a4f0f73SDmitry Karpeev ierr = MPI_Comm_rank(((PetscObject)osm->ois[d])->comm, &srank);CHKERRQ(ierr); 230*c5e4d11fSDmitry Karpeev ierr = PetscViewerGetSubViewer(viewer,((PetscObject)osm->ois[d])->comm, &sviewer);CHKERRQ(ierr); 2316a4f0f73SDmitry Karpeev ierr = ISGetLocalSize(osm->ois[d],&bsz);CHKERRQ(ierr); 2328f3b4b4dSDmitry Karpeev ierr = PetscViewerASCIISynchronizedPrintf(sviewer,"[%d|%d] (subcomm [%d|%d]) local subdomain number %D, local size = %D\n",rank,size,srank,ssize,d,bsz);CHKERRQ(ierr); 2336a4f0f73SDmitry Karpeev ierr = PetscViewerFlush(sviewer);CHKERRQ(ierr); 23406b43e7eSDmitry Karpeev if (view_subdomains) { 23506b43e7eSDmitry Karpeev ierr = PCGASMSubdomainView_Private(pc,d,sviewer);CHKERRQ(ierr); 23606b43e7eSDmitry Karpeev } 2376a4f0f73SDmitry Karpeev if (!pc->setupcalled) { 238*c5e4d11fSDmitry Karpeev ierr = PetscViewerASCIIPrintf(sviewer, "Solver not set up yet: PCSetUp() not yet called\n");CHKERRQ(ierr); 2392fa5cd67SKarl Rupp } else { 24006b43e7eSDmitry Karpeev ierr = KSPView(osm->ksp[d],sviewer);CHKERRQ(ierr); 2416a4f0f73SDmitry Karpeev } 24206b43e7eSDmitry Karpeev ierr = PetscViewerASCIIPrintf(sviewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr); 2437b23a99aSBarry Smith ierr = PetscViewerFlush(sviewer);CHKERRQ(ierr); 244*c5e4d11fSDmitry Karpeev ierr = PetscViewerRestoreSubViewer(viewer,((PetscObject)osm->ois[d])->comm, &sviewer);CHKERRQ(ierr); 24506b43e7eSDmitry Karpeev ++l; 246b1a0a8a3SJed Brown } 24706b43e7eSDmitry Karpeev } 248ce94432eSBarry Smith ierr = MPI_Barrier(PetscObjectComm((PetscObject)pc));CHKERRQ(ierr); 249b1a0a8a3SJed Brown } 2508f3b4b4dSDmitry Karpeev ierr = PetscFree2(numbering,permutation);CHKERRQ(ierr); 251b1a0a8a3SJed Brown ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 252b1a0a8a3SJed Brown } 253*c5e4d11fSDmitry Karpeev ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 254*c5e4d11fSDmitry Karpeev ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 255b1a0a8a3SJed Brown PetscFunctionReturn(0); 256b1a0a8a3SJed Brown } 257b1a0a8a3SJed Brown 2588f3b4b4dSDmitry Karpeev PETSC_INTERN PetscErrorCode PCGASMCreateLocalSubdomains(Mat A, PetscInt nloc, IS *iis[]); 259af538404SDmitry Karpeev 260b1a0a8a3SJed Brown #undef __FUNCT__ 261f746d493SDmitry Karpeev #define __FUNCT__ "PCSetUp_GASM" 262f746d493SDmitry Karpeev static PetscErrorCode PCSetUp_GASM(PC pc) 263b1a0a8a3SJed Brown { 264f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 265b1a0a8a3SJed Brown PetscErrorCode ierr; 266b1a0a8a3SJed Brown PetscBool symset,flg; 26788389755SJed Brown PetscInt i; 268c10bc1a1SDmitry Karpeev PetscMPIInt rank, size; 269b1a0a8a3SJed Brown MatReuse scall = MAT_REUSE_MATRIX; 270b1a0a8a3SJed Brown KSP ksp; 271b1a0a8a3SJed Brown PC subpc; 272b1a0a8a3SJed Brown const char *prefix,*pprefix; 273f746d493SDmitry Karpeev Vec x,y; 2746a4f0f73SDmitry Karpeev PetscInt oni; /* Number of indices in the i-th local outer subdomain. */ 2756a4f0f73SDmitry Karpeev const PetscInt *oidxi; /* Indices from the i-th subdomain local outer subdomain. */ 2766a4f0f73SDmitry Karpeev PetscInt on; /* Number of indices in the disjoint union of local outer subdomains. */ 2776a4f0f73SDmitry Karpeev PetscInt *oidx; /* Indices in the disjoint union of local outer subdomains. */ 2786a4f0f73SDmitry Karpeev IS gois; /* Disjoint union the global indices of outer subdomains. */ 2796a4f0f73SDmitry Karpeev IS goid; /* Identity IS of the size of the disjoint union of outer subdomains. */ 280f746d493SDmitry Karpeev PetscScalar *gxarray, *gyarray; 281930d09c1SFande Kong PetscInt gostart; /* Start of locally-owned indices in the vectors -- osm->gx,osm->gy -- over the disjoint union of outer subdomains. */ 2828f3b4b4dSDmitry Karpeev PetscInt num_subdomains = 0; 2830298fd71SBarry Smith DM *subdomain_dm = NULL; 2848f3b4b4dSDmitry Karpeev char **subdomain_names = NULL; 285f771a274SFande Kong PetscInt *numbering; 2868f3b4b4dSDmitry Karpeev 287b1a0a8a3SJed Brown 288b1a0a8a3SJed Brown PetscFunctionBegin; 289ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr); 290ce94432eSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);CHKERRQ(ierr); 291b1a0a8a3SJed Brown if (!pc->setupcalled) { 292b1a0a8a3SJed Brown 293b1a0a8a3SJed Brown if (!osm->type_set) { 294b1a0a8a3SJed Brown ierr = MatIsSymmetricKnown(pc->pmat,&symset,&flg);CHKERRQ(ierr); 2952fa5cd67SKarl Rupp if (symset && flg) osm->type = PC_GASM_BASIC; 296b1a0a8a3SJed Brown } 297b1a0a8a3SJed Brown 2988f3b4b4dSDmitry Karpeev if (osm->n == PETSC_DETERMINE) { 2998f3b4b4dSDmitry Karpeev if (osm->N != PETSC_DETERMINE) { 3008f3b4b4dSDmitry Karpeev /* No local subdomains given, but the desired number of total subdomains is known, so construct them accordingly. */ 3018f3b4b4dSDmitry Karpeev ierr = PCGASMCreateSubdomains(pc->pmat,osm->N,&osm->n,&osm->iis);CHKERRQ(ierr); 3028f3b4b4dSDmitry Karpeev } else if (osm->dm_subdomains && pc->dm) { 3038f3b4b4dSDmitry Karpeev /* try pc->dm next, if allowed */ 3048f3b4b4dSDmitry Karpeev PetscInt d; 305a35b7b57SDmitry Karpeev IS *inner_subdomain_is, *outer_subdomain_is; 306a35b7b57SDmitry Karpeev ierr = DMCreateDomainDecomposition(pc->dm, &num_subdomains, &subdomain_names, &inner_subdomain_is, &outer_subdomain_is, &subdomain_dm);CHKERRQ(ierr); 307a35b7b57SDmitry Karpeev if (num_subdomains) { 308a35b7b57SDmitry Karpeev ierr = PCGASMSetSubdomains(pc, num_subdomains, inner_subdomain_is, outer_subdomain_is);CHKERRQ(ierr); 30969ca1f37SDmitry Karpeev } 310a35b7b57SDmitry Karpeev for (d = 0; d < num_subdomains; ++d) { 311a35b7b57SDmitry Karpeev if (inner_subdomain_is) {ierr = ISDestroy(&inner_subdomain_is[d]);CHKERRQ(ierr);} 312a35b7b57SDmitry Karpeev if (outer_subdomain_is) {ierr = ISDestroy(&outer_subdomain_is[d]);CHKERRQ(ierr);} 313fdc48646SDmitry Karpeev } 314a35b7b57SDmitry Karpeev ierr = PetscFree(inner_subdomain_is);CHKERRQ(ierr); 315a35b7b57SDmitry Karpeev ierr = PetscFree(outer_subdomain_is);CHKERRQ(ierr); 3168f3b4b4dSDmitry Karpeev } else { 3178f3b4b4dSDmitry Karpeev /* still no subdomains; use one per processor */ 31844fe18e5SDmitry Karpeev osm->nmax = osm->n = 1; 319ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr); 320f746d493SDmitry Karpeev osm->N = size; 3218f3b4b4dSDmitry Karpeev ierr = PCGASMCreateLocalSubdomains(pc->pmat,osm->n,&osm->iis);CHKERRQ(ierr); 322fdc48646SDmitry Karpeev } 32306b43e7eSDmitry Karpeev } 324a35b7b57SDmitry Karpeev if (!osm->iis) { 325a35b7b57SDmitry Karpeev /* 3268f3b4b4dSDmitry Karpeev osm->n was set in PCGASMSetSubdomains(), but the actual subdomains have not been supplied. 3278f3b4b4dSDmitry Karpeev We create the requisite number of local inner subdomains and then expand them into 3288f3b4b4dSDmitry Karpeev out subdomains, if necessary. 329a35b7b57SDmitry Karpeev */ 3308f3b4b4dSDmitry Karpeev ierr = PCGASMCreateLocalSubdomains(pc->pmat,osm->n,&osm->iis);CHKERRQ(ierr); 331a35b7b57SDmitry Karpeev } 3328f3b4b4dSDmitry Karpeev if (!osm->ois) { 3338f3b4b4dSDmitry Karpeev /* 3348f3b4b4dSDmitry Karpeev Initially make outer subdomains the same as inner subdomains. If nonzero additional overlap 3358f3b4b4dSDmitry Karpeev has been requested, copy the inner subdomains over so they can be modified. 3368f3b4b4dSDmitry Karpeev */ 3378f3b4b4dSDmitry Karpeev ierr = PetscMalloc1(osm->n,&osm->ois);CHKERRQ(ierr); 3388f3b4b4dSDmitry Karpeev for (i=0; i<osm->n; ++i) { 3398f3b4b4dSDmitry Karpeev if (osm->overlap > 0) { /* With positive overlap, osm->iis[i] will be modified */ 3408f3b4b4dSDmitry Karpeev ierr = ISDuplicate(osm->iis[i],(osm->ois)+i);CHKERRQ(ierr); 3418f3b4b4dSDmitry Karpeev ierr = ISCopy(osm->iis[i],osm->ois[i]);CHKERRQ(ierr); 3428f3b4b4dSDmitry Karpeev } else { 3438f3b4b4dSDmitry Karpeev ierr = PetscObjectReference((PetscObject)((osm->iis)[i]));CHKERRQ(ierr); 3448f3b4b4dSDmitry Karpeev osm->ois[i] = osm->iis[i]; 3458f3b4b4dSDmitry Karpeev } 3468f3b4b4dSDmitry Karpeev } 3478f3b4b4dSDmitry Karpeev if (osm->overlap > 0) { 3488f3b4b4dSDmitry Karpeev /* Extend the "overlapping" regions by a number of steps */ 349d21f2a47SFande Kong ierr = MatIncreaseOverlapSplit(pc->pmat,osm->n,osm->ois,osm->overlap);CHKERRQ(ierr); 3508f3b4b4dSDmitry Karpeev } 351b1a0a8a3SJed Brown } 3526a4f0f73SDmitry Karpeev 3538f3b4b4dSDmitry Karpeev /* Now the subdomains are defined. Determine their global and max local numbers, if necessary. */ 3548f3b4b4dSDmitry Karpeev if (osm->nmax == PETSC_DETERMINE) { 3558f3b4b4dSDmitry Karpeev PetscMPIInt inwork,outwork; 3568f3b4b4dSDmitry Karpeev /* determine global number of subdomains and the max number of local subdomains */ 3578f3b4b4dSDmitry Karpeev inwork = osm->n; 3588f3b4b4dSDmitry Karpeev ierr = MPI_Allreduce(&inwork,&outwork,1,MPI_INT,MPI_MAX,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr); 3598f3b4b4dSDmitry Karpeev osm->nmax = outwork; 3608f3b4b4dSDmitry Karpeev } 3618f3b4b4dSDmitry Karpeev if (osm->N == PETSC_DETERMINE) { 3628f3b4b4dSDmitry Karpeev /* Determine the number of globally-distinct subdomains and compute a global numbering for them. */ 3638f3b4b4dSDmitry Karpeev ierr = PetscObjectsListGetGlobalNumbering(PetscObjectComm((PetscObject)pc),osm->n,(PetscObject*)osm->ois,&osm->N,NULL);CHKERRQ(ierr); 3648f3b4b4dSDmitry Karpeev } 3658f3b4b4dSDmitry Karpeev 366b1a0a8a3SJed Brown 367b1a0a8a3SJed Brown if (osm->sort_indices) { 368f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 3696a4f0f73SDmitry Karpeev ierr = ISSort(osm->ois[i]);CHKERRQ(ierr); 3706a4f0f73SDmitry Karpeev ierr = ISSort(osm->iis[i]);CHKERRQ(ierr); 371b1a0a8a3SJed Brown } 372b1a0a8a3SJed Brown } 3738f3b4b4dSDmitry Karpeev ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 3748f3b4b4dSDmitry Karpeev ierr = PCGASMPrintSubdomains(pc);CHKERRQ(ierr); 3758f3b4b4dSDmitry Karpeev 3766a4f0f73SDmitry Karpeev /* 3776a4f0f73SDmitry Karpeev Merge the ISs, create merged vectors and restrictions. 3786a4f0f73SDmitry Karpeev */ 3796a4f0f73SDmitry Karpeev /* Merge outer subdomain ISs and construct a restriction onto the disjoint union of local outer subdomains. */ 3806a4f0f73SDmitry Karpeev on = 0; 381f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 3826a4f0f73SDmitry Karpeev ierr = ISGetLocalSize(osm->ois[i],&oni);CHKERRQ(ierr); 3836a4f0f73SDmitry Karpeev on += oni; 384f746d493SDmitry Karpeev } 385785e854fSJed Brown ierr = PetscMalloc1(on, &oidx);CHKERRQ(ierr); 3866a4f0f73SDmitry Karpeev on = 0; 387930d09c1SFande Kong /* Merge local indices together */ 388f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 3896a4f0f73SDmitry Karpeev ierr = ISGetLocalSize(osm->ois[i],&oni);CHKERRQ(ierr); 3906a4f0f73SDmitry Karpeev ierr = ISGetIndices(osm->ois[i],&oidxi);CHKERRQ(ierr); 3916a4f0f73SDmitry Karpeev ierr = PetscMemcpy(oidx+on,oidxi,sizeof(PetscInt)*oni);CHKERRQ(ierr); 3926a4f0f73SDmitry Karpeev ierr = ISRestoreIndices(osm->ois[i],&oidxi);CHKERRQ(ierr); 3936a4f0f73SDmitry Karpeev on += oni; 394f746d493SDmitry Karpeev } 3956a4f0f73SDmitry Karpeev ierr = ISCreateGeneral(((PetscObject)(pc))->comm,on,oidx,PETSC_OWN_POINTER,&gois);CHKERRQ(ierr); 3962a7a6963SBarry Smith ierr = MatCreateVecs(pc->pmat,&x,&y);CHKERRQ(ierr); 397ce94432eSBarry Smith ierr = VecCreateMPI(PetscObjectComm((PetscObject)pc),on,PETSC_DECIDE,&osm->gx);CHKERRQ(ierr); 398f746d493SDmitry Karpeev ierr = VecDuplicate(osm->gx,&osm->gy);CHKERRQ(ierr); 399930d09c1SFande Kong ierr = VecGetOwnershipRange(osm->gx, &gostart, NULL);CHKERRQ(ierr); 400930d09c1SFande Kong ierr = ISCreateStride(PetscObjectComm((PetscObject)pc),on,gostart,1, &goid);CHKERRQ(ierr); 401930d09c1SFande Kong /* gois might indices not on local */ 4026a4f0f73SDmitry Karpeev ierr = VecScatterCreate(x,gois,osm->gx,goid, &(osm->gorestriction));CHKERRQ(ierr); 403f771a274SFande Kong ierr = PetscMalloc1(osm->n,&numbering);CHKERRQ(ierr); 404f771a274SFande Kong ierr = PetscObjectsListGetGlobalNumbering(PetscObjectComm((PetscObject)pc),osm->n,(PetscObject*)osm->ois,NULL,numbering);CHKERRQ(ierr); 4056a4f0f73SDmitry Karpeev ierr = VecDestroy(&x);CHKERRQ(ierr); 4067105d80bSDmitry Karpeev ierr = ISDestroy(&gois);CHKERRQ(ierr); 4072fa5cd67SKarl Rupp 4086a4f0f73SDmitry Karpeev /* Merge inner subdomain ISs and construct a restriction onto the disjoint union of local inner subdomains. */ 4092fa5cd67SKarl Rupp { 4102fa5cd67SKarl Rupp PetscInt ini; /* Number of indices the i-th a local inner subdomain. */ 4116a4f0f73SDmitry Karpeev PetscInt in; /* Number of indices in the disjoint uniont of local inner subdomains. */ 4126a4f0f73SDmitry Karpeev PetscInt *iidx; /* Global indices in the merged local inner subdomain. */ 4136a4f0f73SDmitry Karpeev PetscInt *ioidx; /* Global indices of the disjoint union of inner subdomains within the disjoint union of outer subdomains. */ 4146a4f0f73SDmitry Karpeev IS giis; /* IS for the disjoint union of inner subdomains. */ 4156a4f0f73SDmitry Karpeev IS giois; /* IS for the disjoint union of inner subdomains within the disjoint union of outer subdomains. */ 416f771a274SFande Kong PetscScalar *array; 417f771a274SFande Kong const PetscInt *indices; 418f771a274SFande Kong PetscInt k; 419f746d493SDmitry Karpeev /* */ 4206a4f0f73SDmitry Karpeev on = 0; 421f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 4226a4f0f73SDmitry Karpeev ierr = ISGetLocalSize(osm->ois[i],&oni);CHKERRQ(ierr); 4236a4f0f73SDmitry Karpeev on += oni; 424f746d493SDmitry Karpeev } 425f771a274SFande Kong /* allocate memory */ 426f771a274SFande Kong ierr = PetscMalloc1(on, &iidx);CHKERRQ(ierr); 427f771a274SFande Kong ierr = PetscMalloc1(on, &ioidx);CHKERRQ(ierr); 428f771a274SFande Kong ierr = VecGetArray(y,&array);CHKERRQ(ierr); 429f771a274SFande Kong /* set communicator id */ 430f771a274SFande Kong in = 0; 431f771a274SFande Kong for (i=0; i<osm->n; i++) { 432f771a274SFande Kong ierr = ISGetLocalSize(osm->iis[i],&ini);CHKERRQ(ierr); 433f771a274SFande Kong for (k = 0; k < ini; ++k){ 434f771a274SFande Kong array[in+k] = numbering[i]; 435f771a274SFande Kong } 436f771a274SFande Kong in += ini; 437f771a274SFande Kong } 438f771a274SFande Kong ierr = VecRestoreArray(y,&array);CHKERRQ(ierr); 439f771a274SFande Kong ierr = VecScatterBegin(osm->gorestriction,y,osm->gy,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 440f771a274SFande Kong ierr = VecScatterEnd(osm->gorestriction,y,osm->gy,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 441f771a274SFande Kong ierr = VecGetOwnershipRange(osm->gy,&gostart, NULL);CHKERRQ(ierr); 442f771a274SFande Kong ierr = VecGetArray(osm->gy,&array);CHKERRQ(ierr); 443f771a274SFande Kong on = 0; 444f771a274SFande Kong in = 0; 445f771a274SFande Kong for (i=0; i<osm->n; i++) { 446f771a274SFande Kong ierr = ISGetLocalSize(osm->ois[i],&oni);CHKERRQ(ierr); 447f771a274SFande Kong ierr = ISGetIndices(osm->ois[i],&indices);CHKERRQ(ierr); 448f771a274SFande Kong for (k=0; k<oni; k++) { 449f771a274SFande Kong /* skip overlapping indices */ 45043081b6cSDmitry Karpeev if(PetscRealPart(array[on+k]) != numbering[i]) continue; 451f771a274SFande Kong /* record inner indices */ 452f771a274SFande Kong iidx[in] = indices[k]; 453f771a274SFande Kong ioidx[in++] = gostart+on+k; 454f771a274SFande Kong } 455f771a274SFande Kong ierr = ISRestoreIndices(osm->ois[i], &indices);CHKERRQ(ierr); 456f771a274SFande Kong on += oni; 457f771a274SFande Kong } 458f771a274SFande Kong ierr = VecRestoreArray(osm->gy,&array);CHKERRQ(ierr); 459ce94432eSBarry Smith ierr = ISCreateGeneral(PetscObjectComm((PetscObject)pc),in,iidx,PETSC_OWN_POINTER,&giis);CHKERRQ(ierr); 460ce94432eSBarry Smith ierr = ISCreateGeneral(PetscObjectComm((PetscObject)pc),in,ioidx,PETSC_OWN_POINTER,&giois);CHKERRQ(ierr); 4616a4f0f73SDmitry Karpeev ierr = VecScatterCreate(y,giis,osm->gy,giois,&osm->girestriction);CHKERRQ(ierr); 4626a4f0f73SDmitry Karpeev ierr = VecDestroy(&y);CHKERRQ(ierr); 4636a4f0f73SDmitry Karpeev ierr = ISDestroy(&giis);CHKERRQ(ierr); 4646a4f0f73SDmitry Karpeev ierr = ISDestroy(&giois);CHKERRQ(ierr); 465b1a0a8a3SJed Brown } 4666a4f0f73SDmitry Karpeev ierr = ISDestroy(&goid);CHKERRQ(ierr); 467f771a274SFande Kong ierr = PetscFree(numbering);CHKERRQ(ierr); 4682fa5cd67SKarl Rupp 469f746d493SDmitry Karpeev /* Create the subdomain work vectors. */ 470785e854fSJed Brown ierr = PetscMalloc1(osm->n,&osm->x);CHKERRQ(ierr); 471785e854fSJed Brown ierr = PetscMalloc1(osm->n,&osm->y);CHKERRQ(ierr); 472f746d493SDmitry Karpeev ierr = VecGetArray(osm->gx, &gxarray);CHKERRQ(ierr); 473f746d493SDmitry Karpeev ierr = VecGetArray(osm->gy, &gyarray);CHKERRQ(ierr); 4746a4f0f73SDmitry Karpeev for (i=0, on=0; i<osm->n; ++i, on += oni) { 4756a4f0f73SDmitry Karpeev PetscInt oNi; 4766a4f0f73SDmitry Karpeev ierr = ISGetLocalSize(osm->ois[i],&oni);CHKERRQ(ierr); 477930d09c1SFande Kong /* on a sub communicator */ 4786a4f0f73SDmitry Karpeev ierr = ISGetSize(osm->ois[i],&oNi);CHKERRQ(ierr); 4796a4f0f73SDmitry Karpeev ierr = VecCreateMPIWithArray(((PetscObject)(osm->ois[i]))->comm,1,oni,oNi,gxarray+on,&osm->x[i]);CHKERRQ(ierr); 4806a4f0f73SDmitry Karpeev ierr = VecCreateMPIWithArray(((PetscObject)(osm->ois[i]))->comm,1,oni,oNi,gyarray+on,&osm->y[i]);CHKERRQ(ierr); 481b1a0a8a3SJed Brown } 482f746d493SDmitry Karpeev ierr = VecRestoreArray(osm->gx, &gxarray);CHKERRQ(ierr); 483f746d493SDmitry Karpeev ierr = VecRestoreArray(osm->gy, &gyarray);CHKERRQ(ierr); 4848f3b4b4dSDmitry Karpeev /* Create the subdomain solvers */ 485785e854fSJed Brown ierr = PetscMalloc1(osm->n,&osm->ksp);CHKERRQ(ierr); 4868f3b4b4dSDmitry Karpeev for (i=0; i<osm->n; i++) { 4878f3b4b4dSDmitry Karpeev char subprefix[PETSC_MAX_PATH_LEN+1]; 4886a4f0f73SDmitry Karpeev ierr = KSPCreate(((PetscObject)(osm->ois[i]))->comm,&ksp);CHKERRQ(ierr); 4893bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ksp);CHKERRQ(ierr); 490b1a0a8a3SJed Brown ierr = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,1);CHKERRQ(ierr); 491b1a0a8a3SJed Brown ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr); 4928f3b4b4dSDmitry Karpeev ierr = KSPGetPC(ksp,&subpc);CHKERRQ(ierr); /* Why do we need this here? */ 4938f3b4b4dSDmitry Karpeev if (subdomain_dm) { 4948f3b4b4dSDmitry Karpeev ierr = KSPSetDM(ksp,subdomain_dm[i]);CHKERRQ(ierr); 4958f3b4b4dSDmitry Karpeev ierr = DMDestroy(subdomain_dm+i);CHKERRQ(ierr); 4968f3b4b4dSDmitry Karpeev } 497b1a0a8a3SJed Brown ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 498b1a0a8a3SJed Brown ierr = KSPSetOptionsPrefix(ksp,prefix);CHKERRQ(ierr); 4998f3b4b4dSDmitry Karpeev if (subdomain_names && subdomain_names[i]) { 5008f3b4b4dSDmitry Karpeev ierr = PetscSNPrintf(subprefix,PETSC_MAX_PATH_LEN,"sub_%s_",subdomain_names[i]);CHKERRQ(ierr); 5018f3b4b4dSDmitry Karpeev ierr = KSPAppendOptionsPrefix(ksp,subprefix);CHKERRQ(ierr); 5028f3b4b4dSDmitry Karpeev ierr = PetscFree(subdomain_names[i]);CHKERRQ(ierr); 5038f3b4b4dSDmitry Karpeev } 504b1a0a8a3SJed Brown ierr = KSPAppendOptionsPrefix(ksp,"sub_");CHKERRQ(ierr); 505b1a0a8a3SJed Brown osm->ksp[i] = ksp; 506b1a0a8a3SJed Brown } 5078f3b4b4dSDmitry Karpeev ierr = PetscFree(subdomain_dm);CHKERRQ(ierr); 5088f3b4b4dSDmitry Karpeev ierr = PetscFree(subdomain_names);CHKERRQ(ierr); 509b1a0a8a3SJed Brown scall = MAT_INITIAL_MATRIX; 510b1a0a8a3SJed Brown 5118f3b4b4dSDmitry Karpeev } else { /* if (pc->setupcalled) */ 512b1a0a8a3SJed Brown /* 5138f3b4b4dSDmitry Karpeev Destroy the submatrices from the previous iteration 514b1a0a8a3SJed Brown */ 515b1a0a8a3SJed Brown if (pc->flag == DIFFERENT_NONZERO_PATTERN) { 516f746d493SDmitry Karpeev ierr = MatDestroyMatrices(osm->n,&osm->pmat);CHKERRQ(ierr); 517b1a0a8a3SJed Brown scall = MAT_INITIAL_MATRIX; 518b1a0a8a3SJed Brown } 519b1a0a8a3SJed Brown } 520b1a0a8a3SJed Brown 521b1a0a8a3SJed Brown /* 522f746d493SDmitry Karpeev Extract out the submatrices. 523b1a0a8a3SJed Brown */ 52481a5c4aaSDmitry Karpeev if (size > 1) { 5258f3b4b4dSDmitry Karpeev ierr = MatGetSubMatricesMPI(pc->pmat,osm->n,osm->ois,osm->ois,scall,&osm->pmat);CHKERRQ(ierr); 5262fa5cd67SKarl Rupp } else { 5276a4f0f73SDmitry Karpeev ierr = MatGetSubMatrices(pc->pmat,osm->n,osm->ois,osm->ois,scall,&osm->pmat);CHKERRQ(ierr); 52881a5c4aaSDmitry Karpeev } 529b1a0a8a3SJed Brown if (scall == MAT_INITIAL_MATRIX) { 530b1a0a8a3SJed Brown ierr = PetscObjectGetOptionsPrefix((PetscObject)pc->pmat,&pprefix);CHKERRQ(ierr); 531f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 5323bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)osm->pmat[i]);CHKERRQ(ierr); 533b1a0a8a3SJed Brown ierr = PetscObjectSetOptionsPrefix((PetscObject)osm->pmat[i],pprefix);CHKERRQ(ierr); 534b1a0a8a3SJed Brown } 535b1a0a8a3SJed Brown } 536b1a0a8a3SJed Brown 537b1a0a8a3SJed Brown /* Return control to the user so that the submatrices can be modified (e.g., to apply 538b1a0a8a3SJed Brown different boundary conditions for the submatrices than for the global problem) */ 5396a4f0f73SDmitry Karpeev ierr = PCModifySubMatrices(pc,osm->n,osm->ois,osm->ois,osm->pmat,pc->modifysubmatricesP);CHKERRQ(ierr); 540b1a0a8a3SJed Brown 541b1a0a8a3SJed Brown /* 5426a4f0f73SDmitry Karpeev Loop over submatrices putting them into local ksps 543b1a0a8a3SJed Brown */ 544f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 54523ee1639SBarry Smith ierr = KSPSetOperators(osm->ksp[i],osm->pmat[i],osm->pmat[i]);CHKERRQ(ierr); 546b1a0a8a3SJed Brown if (!pc->setupcalled) { 547b1a0a8a3SJed Brown ierr = KSPSetFromOptions(osm->ksp[i]);CHKERRQ(ierr); 548b1a0a8a3SJed Brown } 549b1a0a8a3SJed Brown } 550b1a0a8a3SJed Brown PetscFunctionReturn(0); 551b1a0a8a3SJed Brown } 552b1a0a8a3SJed Brown 553b1a0a8a3SJed Brown #undef __FUNCT__ 554f746d493SDmitry Karpeev #define __FUNCT__ "PCSetUpOnBlocks_GASM" 555f746d493SDmitry Karpeev static PetscErrorCode PCSetUpOnBlocks_GASM(PC pc) 556b1a0a8a3SJed Brown { 557f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 558b1a0a8a3SJed Brown PetscErrorCode ierr; 559b1a0a8a3SJed Brown PetscInt i; 560b1a0a8a3SJed Brown 561b1a0a8a3SJed Brown PetscFunctionBegin; 562f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 563b1a0a8a3SJed Brown ierr = KSPSetUp(osm->ksp[i]);CHKERRQ(ierr); 564b1a0a8a3SJed Brown } 565b1a0a8a3SJed Brown PetscFunctionReturn(0); 566b1a0a8a3SJed Brown } 567b1a0a8a3SJed Brown 568b1a0a8a3SJed Brown #undef __FUNCT__ 569f746d493SDmitry Karpeev #define __FUNCT__ "PCApply_GASM" 570f746d493SDmitry Karpeev static PetscErrorCode PCApply_GASM(PC pc,Vec x,Vec y) 571b1a0a8a3SJed Brown { 572f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 573b1a0a8a3SJed Brown PetscErrorCode ierr; 574f746d493SDmitry Karpeev PetscInt i; 575b1a0a8a3SJed Brown ScatterMode forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE; 576b1a0a8a3SJed Brown 577b1a0a8a3SJed Brown PetscFunctionBegin; 578b1a0a8a3SJed Brown /* 5796a4f0f73SDmitry Karpeev Support for limiting the restriction or interpolation only to the inner 580b1a0a8a3SJed Brown subdomain values (leaving the other values 0). 581b1a0a8a3SJed Brown */ 582f746d493SDmitry Karpeev if (!(osm->type & PC_GASM_RESTRICT)) { 583b1a0a8a3SJed Brown /* have to zero the work RHS since scatter may leave some slots empty */ 584f746d493SDmitry Karpeev ierr = VecZeroEntries(osm->gx);CHKERRQ(ierr); 5856a4f0f73SDmitry Karpeev ierr = VecScatterBegin(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 5862fa5cd67SKarl Rupp } else { 5876a4f0f73SDmitry Karpeev ierr = VecScatterBegin(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 588b1a0a8a3SJed Brown } 5896a4f0f73SDmitry Karpeev ierr = VecZeroEntries(osm->gy);CHKERRQ(ierr); 5906a4f0f73SDmitry Karpeev if (!(osm->type & PC_GASM_RESTRICT)) { 5916a4f0f73SDmitry Karpeev ierr = VecScatterEnd(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 5922fa5cd67SKarl Rupp } else { 5936a4f0f73SDmitry Karpeev ierr = VecScatterEnd(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 5946a4f0f73SDmitry Karpeev } 59586b96d36SDmitry Karpeev /* do the subdomain solves */ 59686b96d36SDmitry Karpeev for (i=0; i<osm->n; ++i) { 597b1a0a8a3SJed Brown ierr = KSPSolve(osm->ksp[i],osm->x[i],osm->y[i]);CHKERRQ(ierr); 598b1a0a8a3SJed Brown } 599930d09c1SFande Kong /* Do we need to zero y ?? */ 6006a4f0f73SDmitry Karpeev ierr = VecZeroEntries(y);CHKERRQ(ierr); 6016a4f0f73SDmitry Karpeev if (!(osm->type & PC_GASM_INTERPOLATE)) { 6026a4f0f73SDmitry Karpeev ierr = VecScatterBegin(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); 6036a4f0f73SDmitry Karpeev ierr = VecScatterEnd(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); PetscFunctionReturn(0); 6042fa5cd67SKarl Rupp } else { 6056a4f0f73SDmitry Karpeev ierr = VecScatterBegin(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); 6066a4f0f73SDmitry Karpeev ierr = VecScatterEnd(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); PetscFunctionReturn(0); 6076a4f0f73SDmitry Karpeev } 608b1a0a8a3SJed Brown } 609b1a0a8a3SJed Brown 610b1a0a8a3SJed Brown #undef __FUNCT__ 611f746d493SDmitry Karpeev #define __FUNCT__ "PCApplyTranspose_GASM" 612f746d493SDmitry Karpeev static PetscErrorCode PCApplyTranspose_GASM(PC pc,Vec x,Vec y) 613b1a0a8a3SJed Brown { 614f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 615b1a0a8a3SJed Brown PetscErrorCode ierr; 616f746d493SDmitry Karpeev PetscInt i; 617b1a0a8a3SJed Brown ScatterMode forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE; 618b1a0a8a3SJed Brown 619b1a0a8a3SJed Brown PetscFunctionBegin; 620b1a0a8a3SJed Brown /* 621b1a0a8a3SJed Brown Support for limiting the restriction or interpolation to only local 622b1a0a8a3SJed Brown subdomain values (leaving the other values 0). 623b1a0a8a3SJed Brown 624f746d493SDmitry Karpeev Note: these are reversed from the PCApply_GASM() because we are applying the 625b1a0a8a3SJed Brown transpose of the three terms 626b1a0a8a3SJed Brown */ 627f746d493SDmitry Karpeev if (!(osm->type & PC_GASM_INTERPOLATE)) { 628b1a0a8a3SJed Brown /* have to zero the work RHS since scatter may leave some slots empty */ 629f746d493SDmitry Karpeev ierr = VecZeroEntries(osm->gx);CHKERRQ(ierr); 6306a4f0f73SDmitry Karpeev ierr = VecScatterBegin(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 6312fa5cd67SKarl Rupp } else { 6326a4f0f73SDmitry Karpeev ierr = VecScatterBegin(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 633b1a0a8a3SJed Brown } 6346a4f0f73SDmitry Karpeev ierr = VecZeroEntries(osm->gy);CHKERRQ(ierr); 6356a4f0f73SDmitry Karpeev if (!(osm->type & PC_GASM_INTERPOLATE)) { 6366a4f0f73SDmitry Karpeev ierr = VecScatterEnd(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 6372fa5cd67SKarl Rupp } else { 6386a4f0f73SDmitry Karpeev ierr = VecScatterEnd(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 6396a4f0f73SDmitry Karpeev } 640b1a0a8a3SJed Brown /* do the local solves */ 641ab3e923fSDmitry 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. */ 642b1a0a8a3SJed Brown ierr = KSPSolveTranspose(osm->ksp[i],osm->x[i],osm->y[i]);CHKERRQ(ierr); 643b1a0a8a3SJed Brown } 6446a4f0f73SDmitry Karpeev ierr = VecZeroEntries(y);CHKERRQ(ierr); 6456a4f0f73SDmitry Karpeev if (!(osm->type & PC_GASM_RESTRICT)) { 6466a4f0f73SDmitry Karpeev ierr = VecScatterBegin(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); 6476a4f0f73SDmitry Karpeev ierr = VecScatterEnd(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); 6482fa5cd67SKarl Rupp } else { 6496a4f0f73SDmitry Karpeev ierr = VecScatterBegin(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); 6506a4f0f73SDmitry Karpeev ierr = VecScatterEnd(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); 6516a4f0f73SDmitry Karpeev } 652b1a0a8a3SJed Brown PetscFunctionReturn(0); 653b1a0a8a3SJed Brown } 654b1a0a8a3SJed Brown 655b1a0a8a3SJed Brown #undef __FUNCT__ 656a06653b4SBarry Smith #define __FUNCT__ "PCReset_GASM" 657a06653b4SBarry Smith static PetscErrorCode PCReset_GASM(PC pc) 658b1a0a8a3SJed Brown { 659f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 660b1a0a8a3SJed Brown PetscErrorCode ierr; 661b1a0a8a3SJed Brown PetscInt i; 662b1a0a8a3SJed Brown 663b1a0a8a3SJed Brown PetscFunctionBegin; 664b1a0a8a3SJed Brown if (osm->ksp) { 665f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 666a06653b4SBarry Smith ierr = KSPReset(osm->ksp[i]);CHKERRQ(ierr); 667b1a0a8a3SJed Brown } 668b1a0a8a3SJed Brown } 669b1a0a8a3SJed Brown if (osm->pmat) { 670f746d493SDmitry Karpeev if (osm->n > 0) { 671ab3e923fSDmitry Karpeev ierr = MatDestroyMatrices(osm->n,&osm->pmat);CHKERRQ(ierr); 672b1a0a8a3SJed Brown } 673b1a0a8a3SJed Brown } 674d34fcf5fSBarry Smith if (osm->x) { 675f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 676fcfd50ebSBarry Smith ierr = VecDestroy(&osm->x[i]);CHKERRQ(ierr); 677fcfd50ebSBarry Smith ierr = VecDestroy(&osm->y[i]);CHKERRQ(ierr); 678b1a0a8a3SJed Brown } 679d34fcf5fSBarry Smith } 6806bf464f9SBarry Smith ierr = VecDestroy(&osm->gx);CHKERRQ(ierr); 6816bf464f9SBarry Smith ierr = VecDestroy(&osm->gy);CHKERRQ(ierr); 682ab3e923fSDmitry Karpeev 6836a4f0f73SDmitry Karpeev ierr = VecScatterDestroy(&osm->gorestriction);CHKERRQ(ierr); 6846a4f0f73SDmitry Karpeev ierr = VecScatterDestroy(&osm->girestriction);CHKERRQ(ierr); 6858f3b4b4dSDmitry Karpeev if (!osm->user_subdomains) { 6862c112581SDmitry Karpeev ierr = PCGASMDestroySubdomains(osm->n,&osm->ois,&osm->iis);CHKERRQ(ierr); 6878f3b4b4dSDmitry Karpeev osm->N = PETSC_DETERMINE; 6888f3b4b4dSDmitry Karpeev osm->nmax = PETSC_DETERMINE; 6898f3b4b4dSDmitry Karpeev } 690a06653b4SBarry Smith PetscFunctionReturn(0); 691a06653b4SBarry Smith } 692a06653b4SBarry Smith 693a06653b4SBarry Smith #undef __FUNCT__ 694a06653b4SBarry Smith #define __FUNCT__ "PCDestroy_GASM" 695a06653b4SBarry Smith static PetscErrorCode PCDestroy_GASM(PC pc) 696a06653b4SBarry Smith { 697a06653b4SBarry Smith PC_GASM *osm = (PC_GASM*)pc->data; 698a06653b4SBarry Smith PetscErrorCode ierr; 699a06653b4SBarry Smith PetscInt i; 700a06653b4SBarry Smith 701a06653b4SBarry Smith PetscFunctionBegin; 702a06653b4SBarry Smith ierr = PCReset_GASM(pc);CHKERRQ(ierr); 703135757f6SDmitry Karpeev 704135757f6SDmitry Karpeev /* PCReset will not destroy subdomains, if user_subdomains is true. */ 705135757f6SDmitry Karpeev ierr = PCGASMDestroySubdomains(osm->n,&osm->ois,&osm->iis);CHKERRQ(ierr); 706135757f6SDmitry Karpeev 707a06653b4SBarry Smith if (osm->ksp) { 708a06653b4SBarry Smith for (i=0; i<osm->n; i++) { 7096bf464f9SBarry Smith ierr = KSPDestroy(&osm->ksp[i]);CHKERRQ(ierr); 710a06653b4SBarry Smith } 711a06653b4SBarry Smith ierr = PetscFree(osm->ksp);CHKERRQ(ierr); 712a06653b4SBarry Smith } 713a06653b4SBarry Smith ierr = PetscFree(osm->x);CHKERRQ(ierr); 714a06653b4SBarry Smith ierr = PetscFree(osm->y);CHKERRQ(ierr); 715c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 716b1a0a8a3SJed Brown PetscFunctionReturn(0); 717b1a0a8a3SJed Brown } 718b1a0a8a3SJed Brown 719b1a0a8a3SJed Brown #undef __FUNCT__ 720f746d493SDmitry Karpeev #define __FUNCT__ "PCSetFromOptions_GASM" 7218c34d3f5SBarry Smith static PetscErrorCode PCSetFromOptions_GASM(PetscOptions *PetscOptionsObject,PC pc) 722a6dfd86eSKarl Rupp { 723f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 724b1a0a8a3SJed Brown PetscErrorCode ierr; 725b1a0a8a3SJed Brown PetscInt blocks,ovl; 726b1a0a8a3SJed Brown PetscBool symset,flg; 727f746d493SDmitry Karpeev PCGASMType gasmtype; 728b1a0a8a3SJed Brown 729b1a0a8a3SJed Brown PetscFunctionBegin; 730b1a0a8a3SJed Brown /* set the type to symmetric if matrix is symmetric */ 731b1a0a8a3SJed Brown if (!osm->type_set && pc->pmat) { 732b1a0a8a3SJed Brown ierr = MatIsSymmetricKnown(pc->pmat,&symset,&flg);CHKERRQ(ierr); 7332fa5cd67SKarl Rupp if (symset && flg) osm->type = PC_GASM_BASIC; 734b1a0a8a3SJed Brown } 735e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"Generalized additive Schwarz options");CHKERRQ(ierr); 7368f3b4b4dSDmitry 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); 7378f3b4b4dSDmitry Karpeev ierr = PetscOptionsInt("-pc_gasm_total_subdomains","Total number of subdomains across communicator","PCGASMSetTotalSubdomains",osm->N,&blocks,&flg);CHKERRQ(ierr); 73865db9045SDmitry Karpeev if (flg) { 7398f3b4b4dSDmitry Karpeev ierr = PCGASMSetTotalSubdomains(pc,blocks);CHKERRQ(ierr); 74065db9045SDmitry Karpeev } 74106b43e7eSDmitry Karpeev ierr = PetscOptionsInt("-pc_gasm_overlap","Number of overlapping degrees of freedom","PCGASMSetOverlap",osm->overlap,&ovl,&flg);CHKERRQ(ierr); 74265db9045SDmitry Karpeev if (flg) { 74365db9045SDmitry Karpeev ierr = PCGASMSetOverlap(pc,ovl);CHKERRQ(ierr); 744d709ea83SDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 74565db9045SDmitry Karpeev } 746b1a0a8a3SJed Brown flg = PETSC_FALSE; 747f746d493SDmitry Karpeev ierr = PetscOptionsEnum("-pc_gasm_type","Type of restriction/extension","PCGASMSetType",PCGASMTypes,(PetscEnum)osm->type,(PetscEnum*)&gasmtype,&flg);CHKERRQ(ierr); 748f746d493SDmitry Karpeev if (flg) {ierr = PCGASMSetType(pc,gasmtype);CHKERRQ(ierr);} 749b1a0a8a3SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 750b1a0a8a3SJed Brown PetscFunctionReturn(0); 751b1a0a8a3SJed Brown } 752b1a0a8a3SJed Brown 753b1a0a8a3SJed Brown /*------------------------------------------------------------------------------------*/ 754b1a0a8a3SJed Brown 755b1a0a8a3SJed Brown #undef __FUNCT__ 7568f3b4b4dSDmitry Karpeev #define __FUNCT__ "PCGASMSetTotalSubdomains" 7578f3b4b4dSDmitry Karpeev /*@ 7588f3b4b4dSDmitry Karpeev PCGASMSetTotalSubdomains - sets the total number of subdomains to use across the 7598f3b4b4dSDmitry Karpeev communicator. 7608f3b4b4dSDmitry Karpeev Logically collective on pc 7618f3b4b4dSDmitry Karpeev 7628f3b4b4dSDmitry Karpeev Input Parameters: 7638f3b4b4dSDmitry Karpeev + pc - the preconditioner 7648f3b4b4dSDmitry Karpeev - N - total number of subdomains 7658f3b4b4dSDmitry Karpeev 7668f3b4b4dSDmitry Karpeev 7678f3b4b4dSDmitry Karpeev Level: beginner 7688f3b4b4dSDmitry Karpeev 7698f3b4b4dSDmitry Karpeev .keywords: PC, ASM, DM, set, subdomains, additive Schwarz 7708f3b4b4dSDmitry Karpeev 7718f3b4b4dSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMSetOverlap() 7728f3b4b4dSDmitry Karpeev PCGASMCreateSubdomains2D() 7738f3b4b4dSDmitry Karpeev @*/ 7748f3b4b4dSDmitry Karpeev PetscErrorCode PCGASMSetTotalSubdomains(PC pc,PetscInt N) 7758f3b4b4dSDmitry Karpeev { 7768f3b4b4dSDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 7778f3b4b4dSDmitry Karpeev PetscMPIInt size,rank; 7788f3b4b4dSDmitry Karpeev PetscErrorCode ierr; 7798f3b4b4dSDmitry Karpeev 7808f3b4b4dSDmitry Karpeev PetscFunctionBegin; 7818f3b4b4dSDmitry 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); 7828f3b4b4dSDmitry Karpeev if (pc->setupcalled) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetTotalSubdomains() should be called before calling PCSetUp()."); 7838f3b4b4dSDmitry Karpeev 7842c112581SDmitry Karpeev ierr = PCGASMDestroySubdomains(osm->n,&osm->iis,&osm->ois);CHKERRQ(ierr); 7858f3b4b4dSDmitry Karpeev osm->ois = osm->iis = NULL; 7868f3b4b4dSDmitry Karpeev 7878f3b4b4dSDmitry Karpeev ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr); 7888f3b4b4dSDmitry Karpeev ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);CHKERRQ(ierr); 7898f3b4b4dSDmitry Karpeev osm->N = N; 7908f3b4b4dSDmitry Karpeev osm->n = PETSC_DETERMINE; 7918f3b4b4dSDmitry Karpeev osm->nmax = PETSC_DETERMINE; 7928f3b4b4dSDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 7938f3b4b4dSDmitry Karpeev PetscFunctionReturn(0); 7948f3b4b4dSDmitry Karpeev } 7958f3b4b4dSDmitry Karpeev 7968f3b4b4dSDmitry Karpeev #undef __FUNCT__ 79706b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMSetSubdomains_GASM" 798b541e6a4SDmitry Karpeev static PetscErrorCode PCGASMSetSubdomains_GASM(PC pc,PetscInt n,IS iis[],IS ois[]) 799b1a0a8a3SJed Brown { 800f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 801b1a0a8a3SJed Brown PetscErrorCode ierr; 802b1a0a8a3SJed Brown PetscInt i; 803b1a0a8a3SJed Brown 804b1a0a8a3SJed Brown PetscFunctionBegin; 8058f3b4b4dSDmitry Karpeev if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Each process must have 1 or more subdomains, got n = %D",n); 8068f3b4b4dSDmitry Karpeev if (pc->setupcalled) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetSubdomains() should be called before calling PCSetUp()."); 807b1a0a8a3SJed Brown 8082c112581SDmitry Karpeev ierr = PCGASMDestroySubdomains(osm->n,&osm->iis,&osm->ois);CHKERRQ(ierr); 8098f3b4b4dSDmitry Karpeev osm->iis = osm->ois = NULL; 8108f3b4b4dSDmitry Karpeev osm->n = n; 8118f3b4b4dSDmitry Karpeev osm->N = PETSC_DETERMINE; 8128f3b4b4dSDmitry Karpeev osm->nmax = PETSC_DETERMINE; 813a35b7b57SDmitry Karpeev if (ois) { 814785e854fSJed Brown ierr = PetscMalloc1(n,&osm->ois);CHKERRQ(ierr); 8158f3b4b4dSDmitry Karpeev for (i=0; i<n; i++) { 8168f3b4b4dSDmitry Karpeev ierr = PetscObjectReference((PetscObject)ois[i]);CHKERRQ(ierr); 8178f3b4b4dSDmitry Karpeev osm->ois[i] = ois[i]; 8188f3b4b4dSDmitry Karpeev } 8198f3b4b4dSDmitry Karpeev /* 8208f3b4b4dSDmitry Karpeev Since the user set the outer subdomains, even if nontrivial overlap was requested via PCGASMSetOverlap(), 8218f3b4b4dSDmitry Karpeev it will be ignored. To avoid confusion later on (e.g., when viewing the PC), the overlap size is set to -1. 8228f3b4b4dSDmitry Karpeev */ 823b1a0a8a3SJed Brown osm->overlap = -1; 824a35b7b57SDmitry Karpeev if (!iis) { 825785e854fSJed Brown ierr = PetscMalloc1(n,&osm->iis);CHKERRQ(ierr); 826a35b7b57SDmitry Karpeev for (i=0; i<n; i++) { 8278f3b4b4dSDmitry Karpeev for (i=0; i<n; i++) { 8280e8f3646SDmitry Karpeev ierr = PetscObjectReference((PetscObject)ois[i]);CHKERRQ(ierr); 829a35b7b57SDmitry Karpeev osm->iis[i] = ois[i]; 830a35b7b57SDmitry Karpeev } 831a35b7b57SDmitry Karpeev } 832a35b7b57SDmitry Karpeev } 8338f3b4b4dSDmitry Karpeev } 834a35b7b57SDmitry Karpeev if (iis) { 835785e854fSJed Brown ierr = PetscMalloc1(n,&osm->iis);CHKERRQ(ierr); 8368f3b4b4dSDmitry Karpeev for (i=0; i<n; i++) { 8378f3b4b4dSDmitry Karpeev ierr = PetscObjectReference((PetscObject)iis[i]);CHKERRQ(ierr); 8388f3b4b4dSDmitry Karpeev osm->iis[i] = iis[i]; 8398f3b4b4dSDmitry Karpeev } 840a35b7b57SDmitry Karpeev if (!ois) { 841785e854fSJed Brown ierr = PetscMalloc1(n,&osm->ois);CHKERRQ(ierr); 842a35b7b57SDmitry Karpeev for (i=0; i<n; i++) { 8438f3b4b4dSDmitry Karpeev for (i=0; i<n; i++) { 844a35b7b57SDmitry Karpeev ierr = PetscObjectReference((PetscObject)iis[i]);CHKERRQ(ierr); 845a35b7b57SDmitry Karpeev osm->ois[i] = iis[i]; 846a35b7b57SDmitry Karpeev } 8478f3b4b4dSDmitry Karpeev } 8488f3b4b4dSDmitry Karpeev /* If necessary, osm->ois[i] will be expanded using the requested overlap size in PCSetUp_GASM(). */ 849a35b7b57SDmitry Karpeev } 850a35b7b57SDmitry Karpeev } 8518f3b4b4dSDmitry Karpeev if (iis) osm->user_subdomains = PETSC_TRUE; 852b1a0a8a3SJed Brown PetscFunctionReturn(0); 853b1a0a8a3SJed Brown } 854b1a0a8a3SJed Brown 855b1a0a8a3SJed Brown 856b1a0a8a3SJed Brown #undef __FUNCT__ 857f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetOverlap_GASM" 858f7a08781SBarry Smith static PetscErrorCode PCGASMSetOverlap_GASM(PC pc,PetscInt ovl) 859b1a0a8a3SJed Brown { 860f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 861b1a0a8a3SJed Brown 862b1a0a8a3SJed Brown PetscFunctionBegin; 863ce94432eSBarry Smith if (ovl < 0) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap value requested"); 864ce94432eSBarry Smith if (pc->setupcalled && ovl != osm->overlap) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetOverlap() should be called before PCSetUp()."); 8652fa5cd67SKarl Rupp if (!pc->setupcalled) osm->overlap = ovl; 866b1a0a8a3SJed Brown PetscFunctionReturn(0); 867b1a0a8a3SJed Brown } 868b1a0a8a3SJed Brown 869b1a0a8a3SJed Brown #undef __FUNCT__ 870f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetType_GASM" 871f7a08781SBarry Smith static PetscErrorCode PCGASMSetType_GASM(PC pc,PCGASMType type) 872b1a0a8a3SJed Brown { 873f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 874b1a0a8a3SJed Brown 875b1a0a8a3SJed Brown PetscFunctionBegin; 876b1a0a8a3SJed Brown osm->type = type; 877b1a0a8a3SJed Brown osm->type_set = PETSC_TRUE; 878b1a0a8a3SJed Brown PetscFunctionReturn(0); 879b1a0a8a3SJed Brown } 880b1a0a8a3SJed Brown 881b1a0a8a3SJed Brown #undef __FUNCT__ 882f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetSortIndices_GASM" 883f7a08781SBarry Smith static PetscErrorCode PCGASMSetSortIndices_GASM(PC pc,PetscBool doSort) 884b1a0a8a3SJed Brown { 885f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 886b1a0a8a3SJed Brown 887b1a0a8a3SJed Brown PetscFunctionBegin; 888b1a0a8a3SJed Brown osm->sort_indices = doSort; 889b1a0a8a3SJed Brown PetscFunctionReturn(0); 890b1a0a8a3SJed Brown } 891b1a0a8a3SJed Brown 892b1a0a8a3SJed Brown #undef __FUNCT__ 893f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMGetSubKSP_GASM" 89444fe18e5SDmitry Karpeev /* 8958f3b4b4dSDmitry Karpeev FIXME: This routine might need to be modified now that multiple ranks per subdomain are allowed. 89644fe18e5SDmitry Karpeev In particular, it would upset the global subdomain number calculation. 89744fe18e5SDmitry Karpeev */ 898f7a08781SBarry Smith static PetscErrorCode PCGASMGetSubKSP_GASM(PC pc,PetscInt *n,PetscInt *first,KSP **ksp) 899b1a0a8a3SJed Brown { 900f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 901b1a0a8a3SJed Brown PetscErrorCode ierr; 902b1a0a8a3SJed Brown 903b1a0a8a3SJed Brown PetscFunctionBegin; 904ce94432eSBarry 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"); 905b1a0a8a3SJed Brown 9062fa5cd67SKarl Rupp if (n) *n = osm->n; 907ab3e923fSDmitry Karpeev if (first) { 908ce94432eSBarry Smith ierr = MPI_Scan(&osm->n,first,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr); 909ab3e923fSDmitry Karpeev *first -= osm->n; 910b1a0a8a3SJed Brown } 911b1a0a8a3SJed Brown if (ksp) { 912b1a0a8a3SJed Brown /* Assume that local solves are now different; not necessarily 91306b43e7eSDmitry Karpeev true, though! This flag is used only for PCView_GASM() */ 914b1a0a8a3SJed Brown *ksp = osm->ksp; 9156a4f0f73SDmitry Karpeev osm->same_subdomain_solvers = PETSC_FALSE; 916b1a0a8a3SJed Brown } 917b1a0a8a3SJed Brown PetscFunctionReturn(0); 918ab3e923fSDmitry Karpeev } /* PCGASMGetSubKSP_GASM() */ 919b1a0a8a3SJed Brown 920b1a0a8a3SJed Brown #undef __FUNCT__ 92106b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMSetSubdomains" 922b1a0a8a3SJed Brown /*@C 92306b43e7eSDmitry Karpeev PCGASMSetSubdomains - Sets the subdomains for this processor 92406b43e7eSDmitry Karpeev for the additive Schwarz preconditioner. 925b1a0a8a3SJed Brown 9268f3b4b4dSDmitry Karpeev Collective on pc 927b1a0a8a3SJed Brown 928b1a0a8a3SJed Brown Input Parameters: 9298f3b4b4dSDmitry Karpeev + pc - the preconditioner object 93006b43e7eSDmitry Karpeev . n - the number of subdomains for this processor 9318f3b4b4dSDmitry Karpeev . iis - the index sets that define the inner subdomains (or NULL for PETSc to determine subdomains) 9328f3b4b4dSDmitry 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) 933b1a0a8a3SJed Brown 934b1a0a8a3SJed Brown Notes: 93506b43e7eSDmitry Karpeev The IS indices use the parallel, global numbering of the vector entries. 9366a4f0f73SDmitry Karpeev Inner subdomains are those where the correction is applied. 9376a4f0f73SDmitry Karpeev Outer subdomains are those where the residual necessary to obtain the 9386a4f0f73SDmitry Karpeev corrections is obtained (see PCGASMType for the use of inner/outer subdomains). 9396a4f0f73SDmitry Karpeev Both inner and outer subdomains can extend over several processors. 9406a4f0f73SDmitry Karpeev This processor's portion of a subdomain is known as a local subdomain. 9416a4f0f73SDmitry Karpeev 9426a4f0f73SDmitry Karpeev By default the GASM preconditioner uses 1 (local) subdomain per processor. 9436a4f0f73SDmitry Karpeev 944b1a0a8a3SJed Brown 945b1a0a8a3SJed Brown Level: advanced 946b1a0a8a3SJed Brown 94706b43e7eSDmitry Karpeev .keywords: PC, GASM, set, subdomains, additive Schwarz 948b1a0a8a3SJed Brown 9498f3b4b4dSDmitry Karpeev .seealso: PCGASMSetNumSubdomains(), PCGASMSetOverlap(), PCGASMGetSubKSP(), 95006b43e7eSDmitry Karpeev PCGASMCreateSubdomains2D(), PCGASMGetSubdomains() 951b1a0a8a3SJed Brown @*/ 9526a4f0f73SDmitry Karpeev PetscErrorCode PCGASMSetSubdomains(PC pc,PetscInt n,IS iis[],IS ois[]) 953b1a0a8a3SJed Brown { 9548f3b4b4dSDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 955b1a0a8a3SJed Brown PetscErrorCode ierr; 956b1a0a8a3SJed Brown 957b1a0a8a3SJed Brown PetscFunctionBegin; 958b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 9596a4f0f73SDmitry Karpeev ierr = PetscTryMethod(pc,"PCGASMSetSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,n,iis,ois));CHKERRQ(ierr); 9608f3b4b4dSDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 961b1a0a8a3SJed Brown PetscFunctionReturn(0); 962b1a0a8a3SJed Brown } 963b1a0a8a3SJed Brown 964b1a0a8a3SJed Brown 965b1a0a8a3SJed Brown #undef __FUNCT__ 966f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetOverlap" 967b1a0a8a3SJed Brown /*@ 968f746d493SDmitry Karpeev PCGASMSetOverlap - Sets the overlap between a pair of subdomains for the 969b1a0a8a3SJed Brown additive Schwarz preconditioner. Either all or no processors in the 9708f3b4b4dSDmitry Karpeev pc communicator must call this routine. 971b1a0a8a3SJed Brown 9728f3b4b4dSDmitry Karpeev Logically Collective on pc 973b1a0a8a3SJed Brown 974b1a0a8a3SJed Brown Input Parameters: 975b1a0a8a3SJed Brown + pc - the preconditioner context 9768f3b4b4dSDmitry Karpeev - ovl - the amount of overlap between subdomains (ovl >= 0, default value = 0) 977b1a0a8a3SJed Brown 978b1a0a8a3SJed Brown Options Database Key: 97906b43e7eSDmitry Karpeev . -pc_gasm_overlap <overlap> - Sets overlap 980b1a0a8a3SJed Brown 981b1a0a8a3SJed Brown Notes: 98206b43e7eSDmitry Karpeev By default the GASM preconditioner uses 1 subdomain per processor. To use 9838f3b4b4dSDmitry Karpeev multiple subdomain per perocessor or "straddling" subdomains that intersect 9848f3b4b4dSDmitry Karpeev multiple processors use PCGASMSetSubdomains() (or option -pc_gasm_total_subdomains <n>). 985b1a0a8a3SJed Brown 9868f3b4b4dSDmitry Karpeev The overlap defaults to 0, so if one desires that no additional 987b1a0a8a3SJed Brown overlap be computed beyond what may have been set with a call to 9888f3b4b4dSDmitry Karpeev PCGASMSetSubdomains(), then ovl must be set to be 0. In particular, if one does 9898f3b4b4dSDmitry Karpeev not explicitly set the subdomains in application code, then all overlap would be computed 990f746d493SDmitry Karpeev internally by PETSc, and using an overlap of 0 would result in an GASM 991b1a0a8a3SJed Brown variant that is equivalent to the block Jacobi preconditioner. 992b1a0a8a3SJed Brown 993b1a0a8a3SJed Brown Note that one can define initial index sets with any overlap via 99406b43e7eSDmitry Karpeev PCGASMSetSubdomains(); the routine PCGASMSetOverlap() merely allows 99506b43e7eSDmitry Karpeev PETSc to extend that overlap further, if desired. 996b1a0a8a3SJed Brown 997b1a0a8a3SJed Brown Level: intermediate 998b1a0a8a3SJed Brown 999f746d493SDmitry Karpeev .keywords: PC, GASM, set, overlap 1000b1a0a8a3SJed Brown 10018f3b4b4dSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMGetSubKSP(), 100206b43e7eSDmitry Karpeev PCGASMCreateSubdomains2D(), PCGASMGetSubdomains() 1003b1a0a8a3SJed Brown @*/ 10047087cfbeSBarry Smith PetscErrorCode PCGASMSetOverlap(PC pc,PetscInt ovl) 1005b1a0a8a3SJed Brown { 1006b1a0a8a3SJed Brown PetscErrorCode ierr; 10078f3b4b4dSDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 1008b1a0a8a3SJed Brown 1009b1a0a8a3SJed Brown PetscFunctionBegin; 1010b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1011b1a0a8a3SJed Brown PetscValidLogicalCollectiveInt(pc,ovl,2); 1012f746d493SDmitry Karpeev ierr = PetscTryMethod(pc,"PCGASMSetOverlap_C",(PC,PetscInt),(pc,ovl));CHKERRQ(ierr); 10138f3b4b4dSDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 1014b1a0a8a3SJed Brown PetscFunctionReturn(0); 1015b1a0a8a3SJed Brown } 1016b1a0a8a3SJed Brown 1017b1a0a8a3SJed Brown #undef __FUNCT__ 1018f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetType" 1019b1a0a8a3SJed Brown /*@ 1020f746d493SDmitry Karpeev PCGASMSetType - Sets the type of restriction and interpolation used 1021b1a0a8a3SJed Brown for local problems in the additive Schwarz method. 1022b1a0a8a3SJed Brown 1023b1a0a8a3SJed Brown Logically Collective on PC 1024b1a0a8a3SJed Brown 1025b1a0a8a3SJed Brown Input Parameters: 1026b1a0a8a3SJed Brown + pc - the preconditioner context 1027f746d493SDmitry Karpeev - type - variant of GASM, one of 1028b1a0a8a3SJed Brown .vb 1029f746d493SDmitry Karpeev PC_GASM_BASIC - full interpolation and restriction 1030f746d493SDmitry Karpeev PC_GASM_RESTRICT - full restriction, local processor interpolation 1031f746d493SDmitry Karpeev PC_GASM_INTERPOLATE - full interpolation, local processor restriction 1032f746d493SDmitry Karpeev PC_GASM_NONE - local processor restriction and interpolation 1033b1a0a8a3SJed Brown .ve 1034b1a0a8a3SJed Brown 1035b1a0a8a3SJed Brown Options Database Key: 1036f746d493SDmitry Karpeev . -pc_gasm_type [basic,restrict,interpolate,none] - Sets GASM type 1037b1a0a8a3SJed Brown 1038b1a0a8a3SJed Brown Level: intermediate 1039b1a0a8a3SJed Brown 1040f746d493SDmitry Karpeev .keywords: PC, GASM, set, type 1041b1a0a8a3SJed Brown 10428f3b4b4dSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMGetSubKSP(), 1043f746d493SDmitry Karpeev PCGASMCreateSubdomains2D() 1044b1a0a8a3SJed Brown @*/ 10457087cfbeSBarry Smith PetscErrorCode PCGASMSetType(PC pc,PCGASMType type) 1046b1a0a8a3SJed Brown { 1047b1a0a8a3SJed Brown PetscErrorCode ierr; 1048b1a0a8a3SJed Brown 1049b1a0a8a3SJed Brown PetscFunctionBegin; 1050b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1051b1a0a8a3SJed Brown PetscValidLogicalCollectiveEnum(pc,type,2); 1052f746d493SDmitry Karpeev ierr = PetscTryMethod(pc,"PCGASMSetType_C",(PC,PCGASMType),(pc,type));CHKERRQ(ierr); 1053b1a0a8a3SJed Brown PetscFunctionReturn(0); 1054b1a0a8a3SJed Brown } 1055b1a0a8a3SJed Brown 1056b1a0a8a3SJed Brown #undef __FUNCT__ 1057f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetSortIndices" 1058b1a0a8a3SJed Brown /*@ 1059f746d493SDmitry Karpeev PCGASMSetSortIndices - Determines whether subdomain indices are sorted. 1060b1a0a8a3SJed Brown 1061b1a0a8a3SJed Brown Logically Collective on PC 1062b1a0a8a3SJed Brown 1063b1a0a8a3SJed Brown Input Parameters: 1064b1a0a8a3SJed Brown + pc - the preconditioner context 1065b1a0a8a3SJed Brown - doSort - sort the subdomain indices 1066b1a0a8a3SJed Brown 1067b1a0a8a3SJed Brown Level: intermediate 1068b1a0a8a3SJed Brown 1069f746d493SDmitry Karpeev .keywords: PC, GASM, set, type 1070b1a0a8a3SJed Brown 10718f3b4b4dSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMGetSubKSP(), 1072f746d493SDmitry Karpeev PCGASMCreateSubdomains2D() 1073b1a0a8a3SJed Brown @*/ 10747087cfbeSBarry Smith PetscErrorCode PCGASMSetSortIndices(PC pc,PetscBool doSort) 1075b1a0a8a3SJed Brown { 1076b1a0a8a3SJed Brown PetscErrorCode ierr; 1077b1a0a8a3SJed Brown 1078b1a0a8a3SJed Brown PetscFunctionBegin; 1079b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1080b1a0a8a3SJed Brown PetscValidLogicalCollectiveBool(pc,doSort,2); 1081f746d493SDmitry Karpeev ierr = PetscTryMethod(pc,"PCGASMSetSortIndices_C",(PC,PetscBool),(pc,doSort));CHKERRQ(ierr); 1082b1a0a8a3SJed Brown PetscFunctionReturn(0); 1083b1a0a8a3SJed Brown } 1084b1a0a8a3SJed Brown 1085b1a0a8a3SJed Brown #undef __FUNCT__ 1086f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMGetSubKSP" 1087b1a0a8a3SJed Brown /*@C 1088f746d493SDmitry Karpeev PCGASMGetSubKSP - Gets the local KSP contexts for all blocks on 1089b1a0a8a3SJed Brown this processor. 1090b1a0a8a3SJed Brown 1091b1a0a8a3SJed Brown Collective on PC iff first_local is requested 1092b1a0a8a3SJed Brown 1093b1a0a8a3SJed Brown Input Parameter: 1094b1a0a8a3SJed Brown . pc - the preconditioner context 1095b1a0a8a3SJed Brown 1096b1a0a8a3SJed Brown Output Parameters: 10970298fd71SBarry Smith + n_local - the number of blocks on this processor or NULL 10980298fd71SBarry Smith . first_local - the global number of the first block on this processor or NULL, 10990298fd71SBarry Smith all processors must request or all must pass NULL 1100b1a0a8a3SJed Brown - ksp - the array of KSP contexts 1101b1a0a8a3SJed Brown 1102b1a0a8a3SJed Brown Note: 1103f746d493SDmitry Karpeev After PCGASMGetSubKSP() the array of KSPes is not to be freed 1104b1a0a8a3SJed Brown 1105b1a0a8a3SJed Brown Currently for some matrix implementations only 1 block per processor 1106b1a0a8a3SJed Brown is supported. 1107b1a0a8a3SJed Brown 1108f746d493SDmitry Karpeev You must call KSPSetUp() before calling PCGASMGetSubKSP(). 1109b1a0a8a3SJed Brown 1110b1a0a8a3SJed Brown Level: advanced 1111b1a0a8a3SJed Brown 1112f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, get, sub, KSP, context 1113b1a0a8a3SJed Brown 11148f3b4b4dSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMSetOverlap(), 1115f746d493SDmitry Karpeev PCGASMCreateSubdomains2D(), 1116b1a0a8a3SJed Brown @*/ 11177087cfbeSBarry Smith PetscErrorCode PCGASMGetSubKSP(PC pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp[]) 1118b1a0a8a3SJed Brown { 1119b1a0a8a3SJed Brown PetscErrorCode ierr; 1120b1a0a8a3SJed Brown 1121b1a0a8a3SJed Brown PetscFunctionBegin; 1122b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1123f746d493SDmitry Karpeev ierr = PetscUseMethod(pc,"PCGASMGetSubKSP_C",(PC,PetscInt*,PetscInt*,KSP **),(pc,n_local,first_local,ksp));CHKERRQ(ierr); 1124b1a0a8a3SJed Brown PetscFunctionReturn(0); 1125b1a0a8a3SJed Brown } 1126b1a0a8a3SJed Brown 1127b1a0a8a3SJed Brown /* -------------------------------------------------------------------------------------*/ 1128b1a0a8a3SJed Brown /*MC 1129f746d493SDmitry Karpeev PCGASM - Use the (restricted) additive Schwarz method, each block is (approximately) solved with 1130b1a0a8a3SJed Brown its own KSP object. 1131b1a0a8a3SJed Brown 1132b1a0a8a3SJed Brown Options Database Keys: 11338f3b4b4dSDmitry Karpeev + -pc_gasm_total_subdomains <n> - Sets total number of local subdomains to be distributed among processors 113406b43e7eSDmitry Karpeev . -pc_gasm_view_subdomains - activates the printing of subdomain indices in PCView(), -ksp_view or -snes_view 113506b43e7eSDmitry Karpeev . -pc_gasm_print_subdomains - activates the printing of subdomain indices in PCSetUp() 113606b43e7eSDmitry Karpeev . -pc_gasm_overlap <ovl> - Sets overlap by which to (automatically) extend local subdomains 1137f746d493SDmitry Karpeev - -pc_gasm_type [basic,restrict,interpolate,none] - Sets GASM type 1138b1a0a8a3SJed Brown 1139b1a0a8a3SJed Brown IMPORTANT: If you run with, for example, 3 blocks on 1 processor or 3 blocks on 3 processors you 1140f746d493SDmitry Karpeev will get a different convergence rate due to the default option of -pc_gasm_type restrict. Use 1141f746d493SDmitry Karpeev -pc_gasm_type basic to use the standard GASM. 1142b1a0a8a3SJed Brown 1143b1a0a8a3SJed Brown Notes: Each processor can have one or more blocks, but a block cannot be shared by more 1144b1a0a8a3SJed Brown than one processor. Defaults to one block per processor. 1145b1a0a8a3SJed Brown 1146b1a0a8a3SJed Brown To set options on the solvers for each block append -sub_ to all the KSP, and PC 1147b1a0a8a3SJed Brown options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly 1148b1a0a8a3SJed Brown 1149f746d493SDmitry Karpeev To set the options on the solvers separate for each block call PCGASMGetSubKSP() 1150b1a0a8a3SJed Brown and set the options directly on the resulting KSP object (you can access its PC 1151b1a0a8a3SJed Brown with KSPGetPC()) 1152b1a0a8a3SJed Brown 1153b1a0a8a3SJed Brown 1154b1a0a8a3SJed Brown Level: beginner 1155b1a0a8a3SJed Brown 1156b1a0a8a3SJed Brown Concepts: additive Schwarz method 1157b1a0a8a3SJed Brown 1158b1a0a8a3SJed Brown References: 1159b1a0a8a3SJed Brown An additive variant of the Schwarz alternating method for the case of many subregions 1160b1a0a8a3SJed Brown M Dryja, OB Widlund - Courant Institute, New York University Technical report 1161b1a0a8a3SJed Brown 1162b1a0a8a3SJed Brown Domain Decompositions: Parallel Multilevel Methods for Elliptic Partial Differential Equations, 1163b1a0a8a3SJed Brown Barry Smith, Petter Bjorstad, and William Gropp, Cambridge University Press, ISBN 0-521-49589-X. 1164b1a0a8a3SJed Brown 1165b1a0a8a3SJed Brown .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 116649517cdeSBarry Smith PCBJACOBI, PCGASMGetSubKSP(), PCGASMSetSubdomains(), 11678f3b4b4dSDmitry Karpeev PCSetModifySubmatrices(), PCGASMSetOverlap(), PCGASMSetType() 1168b1a0a8a3SJed Brown 1169b1a0a8a3SJed Brown M*/ 1170b1a0a8a3SJed Brown 1171b1a0a8a3SJed Brown #undef __FUNCT__ 1172f746d493SDmitry Karpeev #define __FUNCT__ "PCCreate_GASM" 11738cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_GASM(PC pc) 1174b1a0a8a3SJed Brown { 1175b1a0a8a3SJed Brown PetscErrorCode ierr; 1176f746d493SDmitry Karpeev PC_GASM *osm; 1177b1a0a8a3SJed Brown 1178b1a0a8a3SJed Brown PetscFunctionBegin; 1179b00a9115SJed Brown ierr = PetscNewLog(pc,&osm);CHKERRQ(ierr); 11802fa5cd67SKarl Rupp 11818f3b4b4dSDmitry Karpeev osm->N = PETSC_DETERMINE; 118206b43e7eSDmitry Karpeev osm->n = PETSC_DECIDE; 11838f3b4b4dSDmitry Karpeev osm->nmax = PETSC_DETERMINE; 11848f3b4b4dSDmitry Karpeev osm->overlap = 0; 1185b1a0a8a3SJed Brown osm->ksp = 0; 11866a4f0f73SDmitry Karpeev osm->gorestriction = 0; 11876a4f0f73SDmitry Karpeev osm->girestriction = 0; 1188ab3e923fSDmitry Karpeev osm->gx = 0; 1189ab3e923fSDmitry Karpeev osm->gy = 0; 1190b1a0a8a3SJed Brown osm->x = 0; 1191b1a0a8a3SJed Brown osm->y = 0; 11926a4f0f73SDmitry Karpeev osm->ois = 0; 11936a4f0f73SDmitry Karpeev osm->iis = 0; 1194b1a0a8a3SJed Brown osm->pmat = 0; 1195f746d493SDmitry Karpeev osm->type = PC_GASM_RESTRICT; 11966a4f0f73SDmitry Karpeev osm->same_subdomain_solvers = PETSC_TRUE; 1197b1a0a8a3SJed Brown osm->sort_indices = PETSC_TRUE; 1198d709ea83SDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 1199b1a0a8a3SJed Brown 1200b1a0a8a3SJed Brown pc->data = (void*)osm; 1201f746d493SDmitry Karpeev pc->ops->apply = PCApply_GASM; 1202f746d493SDmitry Karpeev pc->ops->applytranspose = PCApplyTranspose_GASM; 1203f746d493SDmitry Karpeev pc->ops->setup = PCSetUp_GASM; 1204a06653b4SBarry Smith pc->ops->reset = PCReset_GASM; 1205f746d493SDmitry Karpeev pc->ops->destroy = PCDestroy_GASM; 1206f746d493SDmitry Karpeev pc->ops->setfromoptions = PCSetFromOptions_GASM; 1207f746d493SDmitry Karpeev pc->ops->setuponblocks = PCSetUpOnBlocks_GASM; 1208f746d493SDmitry Karpeev pc->ops->view = PCView_GASM; 1209b1a0a8a3SJed Brown pc->ops->applyrichardson = 0; 1210b1a0a8a3SJed Brown 1211bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetSubdomains_C",PCGASMSetSubdomains_GASM);CHKERRQ(ierr); 1212bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetOverlap_C",PCGASMSetOverlap_GASM);CHKERRQ(ierr); 1213bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetType_C",PCGASMSetType_GASM);CHKERRQ(ierr); 1214bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetSortIndices_C",PCGASMSetSortIndices_GASM);CHKERRQ(ierr); 1215bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGASMGetSubKSP_C",PCGASMGetSubKSP_GASM);CHKERRQ(ierr); 1216b1a0a8a3SJed Brown PetscFunctionReturn(0); 1217b1a0a8a3SJed Brown } 1218b1a0a8a3SJed Brown 1219b1a0a8a3SJed Brown 1220b1a0a8a3SJed Brown #undef __FUNCT__ 122106b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMCreateLocalSubdomains" 12228f3b4b4dSDmitry Karpeev PetscErrorCode PCGASMCreateLocalSubdomains(Mat A, PetscInt nloc, IS *iis[]) 1223b1a0a8a3SJed Brown { 1224b1a0a8a3SJed Brown MatPartitioning mpart; 1225b1a0a8a3SJed Brown const char *prefix; 122611bd1e4dSLisandro Dalcin PetscErrorCode (*f)(Mat,MatReuse,Mat*); 1227b1a0a8a3SJed Brown PetscMPIInt size; 1228b1a0a8a3SJed Brown PetscInt i,j,rstart,rend,bs; 122911bd1e4dSLisandro Dalcin PetscBool isbaij = PETSC_FALSE,foundpart = PETSC_FALSE; 12300298fd71SBarry Smith Mat Ad = NULL, adj; 1231b1a0a8a3SJed Brown IS ispart,isnumb,*is; 1232b1a0a8a3SJed Brown PetscErrorCode ierr; 1233b1a0a8a3SJed Brown 1234b1a0a8a3SJed Brown PetscFunctionBegin; 12358f3b4b4dSDmitry Karpeev if (nloc < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of local subdomains must > 0, got nloc = %D",nloc); 1236b1a0a8a3SJed Brown 1237b1a0a8a3SJed Brown /* Get prefix, row distribution, and block size */ 1238b1a0a8a3SJed Brown ierr = MatGetOptionsPrefix(A,&prefix);CHKERRQ(ierr); 1239b1a0a8a3SJed Brown ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 1240b1a0a8a3SJed Brown ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr); 1241b1a0a8a3SJed 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); 1242b1a0a8a3SJed Brown 1243b1a0a8a3SJed Brown /* Get diagonal block from matrix if possible */ 1244ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr); 12450005d66cSJed Brown ierr = PetscObjectQueryFunction((PetscObject)A,"MatGetDiagonalBlock_C",&f);CHKERRQ(ierr); 1246b1a0a8a3SJed Brown if (f) { 124711bd1e4dSLisandro Dalcin ierr = MatGetDiagonalBlock(A,&Ad);CHKERRQ(ierr); 1248b1a0a8a3SJed Brown } else if (size == 1) { 124911bd1e4dSLisandro Dalcin Ad = A; 1250b1a0a8a3SJed Brown } 1251b1a0a8a3SJed Brown if (Ad) { 1252251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)Ad,MATSEQBAIJ,&isbaij);CHKERRQ(ierr); 1253251f4c67SDmitry Karpeev if (!isbaij) {ierr = PetscObjectTypeCompare((PetscObject)Ad,MATSEQSBAIJ,&isbaij);CHKERRQ(ierr);} 1254b1a0a8a3SJed Brown } 12558f3b4b4dSDmitry Karpeev if (Ad && nloc > 1) { 1256b1a0a8a3SJed Brown PetscBool match,done; 1257b1a0a8a3SJed Brown /* Try to setup a good matrix partitioning if available */ 1258b1a0a8a3SJed Brown ierr = MatPartitioningCreate(PETSC_COMM_SELF,&mpart);CHKERRQ(ierr); 1259b1a0a8a3SJed Brown ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr); 1260b1a0a8a3SJed Brown ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr); 1261251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGCURRENT,&match);CHKERRQ(ierr); 1262b1a0a8a3SJed Brown if (!match) { 1263251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGSQUARE,&match);CHKERRQ(ierr); 1264b1a0a8a3SJed Brown } 1265b1a0a8a3SJed Brown if (!match) { /* assume a "good" partitioner is available */ 12661a83f524SJed Brown PetscInt na; 12671a83f524SJed Brown const PetscInt *ia,*ja; 1268b1a0a8a3SJed Brown ierr = MatGetRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr); 1269b1a0a8a3SJed Brown if (done) { 1270b1a0a8a3SJed Brown /* Build adjacency matrix by hand. Unfortunately a call to 1271b1a0a8a3SJed Brown MatConvert(Ad,MATMPIADJ,MAT_INITIAL_MATRIX,&adj) will 1272b1a0a8a3SJed Brown remove the block-aij structure and we cannot expect 1273b1a0a8a3SJed Brown MatPartitioning to split vertices as we need */ 12741a83f524SJed Brown PetscInt i,j,len,nnz,cnt,*iia=0,*jja=0; 12751a83f524SJed Brown const PetscInt *row; 1276b1a0a8a3SJed Brown nnz = 0; 1277b1a0a8a3SJed Brown for (i=0; i<na; i++) { /* count number of nonzeros */ 1278b1a0a8a3SJed Brown len = ia[i+1] - ia[i]; 1279b1a0a8a3SJed Brown row = ja + ia[i]; 1280b1a0a8a3SJed Brown for (j=0; j<len; j++) { 1281b1a0a8a3SJed Brown if (row[j] == i) { /* don't count diagonal */ 1282b1a0a8a3SJed Brown len--; break; 1283b1a0a8a3SJed Brown } 1284b1a0a8a3SJed Brown } 1285b1a0a8a3SJed Brown nnz += len; 1286b1a0a8a3SJed Brown } 1287854ce69bSBarry Smith ierr = PetscMalloc1(na+1,&iia);CHKERRQ(ierr); 1288854ce69bSBarry Smith ierr = PetscMalloc1(nnz,&jja);CHKERRQ(ierr); 1289b1a0a8a3SJed Brown nnz = 0; 1290b1a0a8a3SJed Brown iia[0] = 0; 1291b1a0a8a3SJed Brown for (i=0; i<na; i++) { /* fill adjacency */ 1292b1a0a8a3SJed Brown cnt = 0; 1293b1a0a8a3SJed Brown len = ia[i+1] - ia[i]; 1294b1a0a8a3SJed Brown row = ja + ia[i]; 1295b1a0a8a3SJed Brown for (j=0; j<len; j++) { 12962fa5cd67SKarl Rupp if (row[j] != i) jja[nnz+cnt++] = row[j]; /* if not diagonal */ 1297b1a0a8a3SJed Brown } 1298b1a0a8a3SJed Brown nnz += cnt; 1299b1a0a8a3SJed Brown iia[i+1] = nnz; 1300b1a0a8a3SJed Brown } 1301b1a0a8a3SJed Brown /* Partitioning of the adjacency matrix */ 13020298fd71SBarry Smith ierr = MatCreateMPIAdj(PETSC_COMM_SELF,na,na,iia,jja,NULL,&adj);CHKERRQ(ierr); 1303b1a0a8a3SJed Brown ierr = MatPartitioningSetAdjacency(mpart,adj);CHKERRQ(ierr); 13048f3b4b4dSDmitry Karpeev ierr = MatPartitioningSetNParts(mpart,nloc);CHKERRQ(ierr); 1305b1a0a8a3SJed Brown ierr = MatPartitioningApply(mpart,&ispart);CHKERRQ(ierr); 1306b1a0a8a3SJed Brown ierr = ISPartitioningToNumbering(ispart,&isnumb);CHKERRQ(ierr); 13076bf464f9SBarry Smith ierr = MatDestroy(&adj);CHKERRQ(ierr); 1308b1a0a8a3SJed Brown foundpart = PETSC_TRUE; 1309b1a0a8a3SJed Brown } 1310b1a0a8a3SJed Brown ierr = MatRestoreRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr); 1311b1a0a8a3SJed Brown } 13126bf464f9SBarry Smith ierr = MatPartitioningDestroy(&mpart);CHKERRQ(ierr); 1313b1a0a8a3SJed Brown } 13148f3b4b4dSDmitry Karpeev ierr = PetscMalloc1(nloc,&is);CHKERRQ(ierr); 1315b1a0a8a3SJed Brown if (!foundpart) { 1316b1a0a8a3SJed Brown 1317b1a0a8a3SJed Brown /* Partitioning by contiguous chunks of rows */ 1318b1a0a8a3SJed Brown 1319b1a0a8a3SJed Brown PetscInt mbs = (rend-rstart)/bs; 1320b1a0a8a3SJed Brown PetscInt start = rstart; 13218f3b4b4dSDmitry Karpeev for (i=0; i<nloc; i++) { 13228f3b4b4dSDmitry Karpeev PetscInt count = (mbs/nloc + ((mbs % nloc) > i)) * bs; 1323b1a0a8a3SJed Brown ierr = ISCreateStride(PETSC_COMM_SELF,count,start,1,&is[i]);CHKERRQ(ierr); 1324b1a0a8a3SJed Brown start += count; 1325b1a0a8a3SJed Brown } 1326b1a0a8a3SJed Brown 1327b1a0a8a3SJed Brown } else { 1328b1a0a8a3SJed Brown 1329b1a0a8a3SJed Brown /* Partitioning by adjacency of diagonal block */ 1330b1a0a8a3SJed Brown 1331b1a0a8a3SJed Brown const PetscInt *numbering; 1332b1a0a8a3SJed Brown PetscInt *count,nidx,*indices,*newidx,start=0; 1333b1a0a8a3SJed Brown /* Get node count in each partition */ 13348f3b4b4dSDmitry Karpeev ierr = PetscMalloc1(nloc,&count);CHKERRQ(ierr); 13358f3b4b4dSDmitry Karpeev ierr = ISPartitioningCount(ispart,nloc,count);CHKERRQ(ierr); 1336b1a0a8a3SJed Brown if (isbaij && bs > 1) { /* adjust for the block-aij case */ 13378f3b4b4dSDmitry Karpeev for (i=0; i<nloc; i++) count[i] *= bs; 1338b1a0a8a3SJed Brown } 1339b1a0a8a3SJed Brown /* Build indices from node numbering */ 1340b1a0a8a3SJed Brown ierr = ISGetLocalSize(isnumb,&nidx);CHKERRQ(ierr); 1341785e854fSJed Brown ierr = PetscMalloc1(nidx,&indices);CHKERRQ(ierr); 1342b1a0a8a3SJed Brown for (i=0; i<nidx; i++) indices[i] = i; /* needs to be initialized */ 1343b1a0a8a3SJed Brown ierr = ISGetIndices(isnumb,&numbering);CHKERRQ(ierr); 1344b1a0a8a3SJed Brown ierr = PetscSortIntWithPermutation(nidx,numbering,indices);CHKERRQ(ierr); 1345b1a0a8a3SJed Brown ierr = ISRestoreIndices(isnumb,&numbering);CHKERRQ(ierr); 1346b1a0a8a3SJed Brown if (isbaij && bs > 1) { /* adjust for the block-aij case */ 1347785e854fSJed Brown ierr = PetscMalloc1(nidx*bs,&newidx);CHKERRQ(ierr); 13482fa5cd67SKarl Rupp for (i=0; i<nidx; i++) { 13492fa5cd67SKarl Rupp for (j=0; j<bs; j++) newidx[i*bs+j] = indices[i]*bs + j; 13502fa5cd67SKarl Rupp } 1351b1a0a8a3SJed Brown ierr = PetscFree(indices);CHKERRQ(ierr); 1352b1a0a8a3SJed Brown nidx *= bs; 1353b1a0a8a3SJed Brown indices = newidx; 1354b1a0a8a3SJed Brown } 1355b1a0a8a3SJed Brown /* Shift to get global indices */ 1356b1a0a8a3SJed Brown for (i=0; i<nidx; i++) indices[i] += rstart; 1357b1a0a8a3SJed Brown 1358b1a0a8a3SJed Brown /* Build the index sets for each block */ 13598f3b4b4dSDmitry Karpeev for (i=0; i<nloc; i++) { 1360b1a0a8a3SJed Brown ierr = ISCreateGeneral(PETSC_COMM_SELF,count[i],&indices[start],PETSC_COPY_VALUES,&is[i]);CHKERRQ(ierr); 1361b1a0a8a3SJed Brown ierr = ISSort(is[i]);CHKERRQ(ierr); 1362b1a0a8a3SJed Brown start += count[i]; 1363b1a0a8a3SJed Brown } 1364b1a0a8a3SJed Brown 1365302440fdSBarry Smith ierr = PetscFree(count);CHKERRQ(ierr); 1366302440fdSBarry Smith ierr = PetscFree(indices);CHKERRQ(ierr); 1367fcfd50ebSBarry Smith ierr = ISDestroy(&isnumb);CHKERRQ(ierr); 1368fcfd50ebSBarry Smith ierr = ISDestroy(&ispart);CHKERRQ(ierr); 1369b1a0a8a3SJed Brown } 13706a4f0f73SDmitry Karpeev *iis = is; 13718f3b4b4dSDmitry Karpeev PetscFunctionReturn(0); 13728f3b4b4dSDmitry Karpeev } 13738f3b4b4dSDmitry Karpeev 13748f3b4b4dSDmitry Karpeev #undef __FUNCT__ 13758f3b4b4dSDmitry Karpeev #define __FUNCT__ "PCGASMCreateStraddlingSubdomains" 1376b541e6a4SDmitry Karpeev PETSC_INTERN PetscErrorCode PCGASMCreateStraddlingSubdomains(Mat A,PetscInt N,PetscInt *n,IS *iis[]) 13778f3b4b4dSDmitry Karpeev { 13788f3b4b4dSDmitry Karpeev PetscErrorCode ierr; 13798f3b4b4dSDmitry Karpeev 13808f3b4b4dSDmitry Karpeev PetscFunctionBegin; 13818f3b4b4dSDmitry Karpeev ierr = MatSubdomainsCreateCoalesce(A,N,n,iis);CHKERRQ(ierr); 13828f3b4b4dSDmitry Karpeev PetscFunctionReturn(0); 13838f3b4b4dSDmitry Karpeev } 13848f3b4b4dSDmitry Karpeev 13858f3b4b4dSDmitry Karpeev 13868f3b4b4dSDmitry Karpeev 13878f3b4b4dSDmitry Karpeev #undef __FUNCT__ 13888f3b4b4dSDmitry Karpeev #define __FUNCT__ "PCGASMCreateSubdomains" 13898f3b4b4dSDmitry Karpeev /*@C 13908f3b4b4dSDmitry Karpeev PCGASMCreateSubdomains - Creates n index sets defining n nonoverlapping subdomains for the additive 13918f3b4b4dSDmitry Karpeev Schwarz preconditioner for a any problem based on its matrix. 13928f3b4b4dSDmitry Karpeev 13938f3b4b4dSDmitry Karpeev Collective 13948f3b4b4dSDmitry Karpeev 13958f3b4b4dSDmitry Karpeev Input Parameters: 13968f3b4b4dSDmitry Karpeev + A - The global matrix operator 13978f3b4b4dSDmitry Karpeev - N - the number of global subdomains requested 13988f3b4b4dSDmitry Karpeev 13998f3b4b4dSDmitry Karpeev Output Parameters: 14008f3b4b4dSDmitry Karpeev + n - the number of subdomains created on this processor 14018f3b4b4dSDmitry Karpeev - iis - the array of index sets defining the local inner subdomains (on which the correction is applied) 14028f3b4b4dSDmitry Karpeev 14038f3b4b4dSDmitry Karpeev Level: advanced 14048f3b4b4dSDmitry Karpeev 14058f3b4b4dSDmitry Karpeev Note: When N >= A's communicator size, each subdomain is local -- contained within a single processor. 14068f3b4b4dSDmitry Karpeev When N < size, the subdomains are 'straddling' (processor boundaries) and are no longer local. 14078f3b4b4dSDmitry Karpeev The resulting subdomains can be use in PCGASMSetSubdomains(pc,n,iss,NULL). The overlapping 14088f3b4b4dSDmitry Karpeev outer subdomains will be automatically generated from these according to the requested amount of 14098f3b4b4dSDmitry Karpeev overlap; this is currently supported only with local subdomains. 14108f3b4b4dSDmitry Karpeev 14118f3b4b4dSDmitry Karpeev 14128f3b4b4dSDmitry Karpeev .keywords: PC, GASM, additive Schwarz, create, subdomains, unstructured grid 14138f3b4b4dSDmitry Karpeev 14148f3b4b4dSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMDestroySubdomains() 14158f3b4b4dSDmitry Karpeev @*/ 1416b541e6a4SDmitry Karpeev PetscErrorCode PCGASMCreateSubdomains(Mat A,PetscInt N,PetscInt *n,IS *iis[]) 14178f3b4b4dSDmitry Karpeev { 14188f3b4b4dSDmitry Karpeev PetscMPIInt size; 14198f3b4b4dSDmitry Karpeev PetscErrorCode ierr; 14208f3b4b4dSDmitry Karpeev 14218f3b4b4dSDmitry Karpeev PetscFunctionBegin; 14228f3b4b4dSDmitry Karpeev PetscValidHeaderSpecific(A,MAT_CLASSID,1); 14238f3b4b4dSDmitry Karpeev PetscValidPointer(iis,4); 14248f3b4b4dSDmitry Karpeev 14258f3b4b4dSDmitry Karpeev if (N < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Number of subdomains must be > 0, N = %D",N); 14268f3b4b4dSDmitry Karpeev ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr); 14278f3b4b4dSDmitry Karpeev if (N >= size) { 14288f3b4b4dSDmitry Karpeev *n = N/size + (N%size); 14298f3b4b4dSDmitry Karpeev ierr = PCGASMCreateLocalSubdomains(A,*n,iis);CHKERRQ(ierr); 14306a4f0f73SDmitry Karpeev } else { 14318f3b4b4dSDmitry Karpeev ierr = PCGASMCreateStraddlingSubdomains(A,N,n,iis);CHKERRQ(ierr); 14326a4f0f73SDmitry Karpeev } 1433b1a0a8a3SJed Brown PetscFunctionReturn(0); 1434b1a0a8a3SJed Brown } 1435b1a0a8a3SJed Brown 1436b1a0a8a3SJed Brown #undef __FUNCT__ 1437f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMDestroySubdomains" 1438b1a0a8a3SJed Brown /*@C 1439f746d493SDmitry Karpeev PCGASMDestroySubdomains - Destroys the index sets created with 14408f3b4b4dSDmitry Karpeev PCGASMCreateSubdomains() or PCGASMCreateSubdomains2D. Should be 144106b43e7eSDmitry Karpeev called after setting subdomains with PCGASMSetSubdomains(). 1442b1a0a8a3SJed Brown 1443b1a0a8a3SJed Brown Collective 1444b1a0a8a3SJed Brown 1445b1a0a8a3SJed Brown Input Parameters: 1446b1a0a8a3SJed Brown + n - the number of index sets 14476a4f0f73SDmitry Karpeev . iis - the array of inner subdomains, 14480298fd71SBarry Smith - ois - the array of outer subdomains, can be NULL 1449b1a0a8a3SJed Brown 14506a4f0f73SDmitry Karpeev Level: intermediate 14516a4f0f73SDmitry Karpeev 14526a4f0f73SDmitry Karpeev Notes: this is merely a convenience subroutine that walks each list, 14532c112581SDmitry Karpeev destroys each IS on the list, and then frees the list. At the end the 14542c112581SDmitry Karpeev list pointers are set to NULL. 1455b1a0a8a3SJed Brown 1456f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, create, subdomains, unstructured grid 1457b1a0a8a3SJed Brown 14588f3b4b4dSDmitry Karpeev .seealso: PCGASMCreateSubdomains(), PCGASMSetSubdomains() 1459b1a0a8a3SJed Brown @*/ 14602c112581SDmitry Karpeev PetscErrorCode PCGASMDestroySubdomains(PetscInt n,IS **iis,IS **ois) 1461b1a0a8a3SJed Brown { 1462b1a0a8a3SJed Brown PetscInt i; 1463b1a0a8a3SJed Brown PetscErrorCode ierr; 14645fd66863SKarl Rupp 1465b1a0a8a3SJed Brown PetscFunctionBegin; 1466a35b7b57SDmitry Karpeev if (n <= 0) PetscFunctionReturn(0); 14676a4f0f73SDmitry Karpeev if (ois) { 14682c112581SDmitry Karpeev PetscValidPointer(ois,3); 14692c112581SDmitry Karpeev if (*ois) { 14702c112581SDmitry Karpeev PetscValidPointer(*ois,3); 1471a35b7b57SDmitry Karpeev for (i=0; i<n; i++) { 14722c112581SDmitry Karpeev ierr = ISDestroy(&(*ois)[i]);CHKERRQ(ierr); 1473a35b7b57SDmitry Karpeev } 1474135757f6SDmitry Karpeev ierr = PetscFree((*ois));CHKERRQ(ierr); 14752c112581SDmitry Karpeev } 1476b1a0a8a3SJed Brown } 1477b9d0fdaaSFande Kong if (iis) { 1478b9d0fdaaSFande Kong PetscValidPointer(iis,2); 1479b9d0fdaaSFande Kong if (*iis) { 1480b9d0fdaaSFande Kong PetscValidPointer(*iis,2); 1481b9d0fdaaSFande Kong for (i=0; i<n; i++) { 1482b9d0fdaaSFande Kong ierr = ISDestroy(&(*iis)[i]);CHKERRQ(ierr); 1483b9d0fdaaSFande Kong } 1484b9d0fdaaSFande Kong ierr = PetscFree((*iis));CHKERRQ(ierr); 1485b9d0fdaaSFande Kong } 1486b9d0fdaaSFande Kong } 1487b1a0a8a3SJed Brown PetscFunctionReturn(0); 1488b1a0a8a3SJed Brown } 1489b1a0a8a3SJed Brown 1490af538404SDmitry Karpeev 1491af538404SDmitry Karpeev #define PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,xleft_loc,ylow_loc,xright_loc,yhigh_loc,n) \ 1492af538404SDmitry Karpeev { \ 1493af538404SDmitry Karpeev PetscInt first_row = first/M, last_row = last/M+1; \ 1494af538404SDmitry Karpeev /* \ 1495af538404SDmitry Karpeev Compute ylow_loc and yhigh_loc so that (ylow_loc,xleft) and (yhigh_loc,xright) are the corners \ 1496af538404SDmitry Karpeev of the bounding box of the intersection of the subdomain with the local ownership range (local \ 1497af538404SDmitry Karpeev subdomain). \ 1498af538404SDmitry Karpeev Also compute xleft_loc and xright_loc as the lower and upper bounds on the first and last rows \ 1499af538404SDmitry Karpeev of the intersection. \ 1500af538404SDmitry Karpeev */ \ 1501af538404SDmitry Karpeev /* ylow_loc is the grid row containing the first element of the local sumbdomain */ \ 1502eec7e3faSDmitry Karpeev *ylow_loc = PetscMax(first_row,ylow); \ 1503af538404SDmitry Karpeev /* xleft_loc is the offset of first element of the local subdomain within its grid row (might actually be outside the local subdomain) */ \ 1504eec7e3faSDmitry Karpeev *xleft_loc = *ylow_loc==first_row ? PetscMax(first%M,xleft) : xleft; \ 1505af538404SDmitry Karpeev /* yhigh_loc is the grid row above the last local subdomain element */ \ 1506eec7e3faSDmitry Karpeev *yhigh_loc = PetscMin(last_row,yhigh); \ 1507af538404SDmitry Karpeev /* xright is the offset of the end of the local subdomain within its grid row (might actually be outside the local subdomain) */ \ 1508eec7e3faSDmitry Karpeev *xright_loc = *yhigh_loc==last_row ? PetscMin(xright,last%M) : xright; \ 1509af538404SDmitry Karpeev /* Now compute the size of the local subdomain n. */ \ 1510c3518bceSDmitry Karpeev *n = 0; \ 1511eec7e3faSDmitry Karpeev if (*ylow_loc < *yhigh_loc) { \ 1512af538404SDmitry Karpeev PetscInt width = xright-xleft; \ 1513eec7e3faSDmitry Karpeev *n += width*(*yhigh_loc-*ylow_loc-1); \ 1514eec7e3faSDmitry Karpeev *n += PetscMin(PetscMax(*xright_loc-xleft,0),width); \ 1515eec7e3faSDmitry Karpeev *n -= PetscMin(PetscMax(*xleft_loc-xleft,0), width); \ 1516af538404SDmitry Karpeev } \ 1517af538404SDmitry Karpeev } 1518af538404SDmitry Karpeev 1519c3518bceSDmitry Karpeev 1520eec7e3faSDmitry Karpeev 1521eec7e3faSDmitry Karpeev #undef __FUNCT__ 1522f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMCreateSubdomains2D" 1523b1a0a8a3SJed Brown /*@ 1524f746d493SDmitry Karpeev PCGASMCreateSubdomains2D - Creates the index sets for the overlapping Schwarz 1525b1a0a8a3SJed Brown preconditioner for a two-dimensional problem on a regular grid. 1526b1a0a8a3SJed Brown 1527af538404SDmitry Karpeev Collective 1528b1a0a8a3SJed Brown 1529b1a0a8a3SJed Brown Input Parameters: 153006b43e7eSDmitry Karpeev + M, N - the global number of grid points in the x and y directions 1531af538404SDmitry Karpeev . Mdomains, Ndomains - the global number of subdomains in the x and y directions 1532b1a0a8a3SJed Brown . dof - degrees of freedom per node 1533b1a0a8a3SJed Brown - overlap - overlap in mesh lines 1534b1a0a8a3SJed Brown 1535b1a0a8a3SJed Brown Output Parameters: 1536af538404SDmitry Karpeev + Nsub - the number of local subdomains created 15376a4f0f73SDmitry Karpeev . iis - array of index sets defining inner (nonoverlapping) subdomains 15386a4f0f73SDmitry Karpeev - ois - array of index sets defining outer (overlapping, if overlap > 0) subdomains 1539b1a0a8a3SJed Brown 1540b1a0a8a3SJed Brown 1541b1a0a8a3SJed Brown Level: advanced 1542b1a0a8a3SJed Brown 1543f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, create, subdomains, 2D, regular grid 1544b1a0a8a3SJed Brown 15458f3b4b4dSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMGetSubKSP(), PCGASMSetOverlap() 1546b1a0a8a3SJed Brown @*/ 15476a4f0f73SDmitry Karpeev PetscErrorCode PCGASMCreateSubdomains2D(PC pc,PetscInt M,PetscInt N,PetscInt Mdomains,PetscInt Ndomains,PetscInt dof,PetscInt overlap,PetscInt *nsub,IS **iis,IS **ois) 1548b1a0a8a3SJed Brown { 1549b1a0a8a3SJed Brown PetscErrorCode ierr; 15502bbb417fSDmitry Karpeev PetscMPIInt size, rank; 15512bbb417fSDmitry Karpeev PetscInt i, j; 15522bbb417fSDmitry Karpeev PetscInt maxheight, maxwidth; 15532bbb417fSDmitry Karpeev PetscInt xstart, xleft, xright, xleft_loc, xright_loc; 15542bbb417fSDmitry Karpeev PetscInt ystart, ylow, yhigh, ylow_loc, yhigh_loc; 1555eec7e3faSDmitry Karpeev PetscInt x[2][2], y[2][2], n[2]; 15562bbb417fSDmitry Karpeev PetscInt first, last; 15572bbb417fSDmitry Karpeev PetscInt nidx, *idx; 15582bbb417fSDmitry Karpeev PetscInt ii,jj,s,q,d; 1559af538404SDmitry Karpeev PetscInt k,kk; 15602bbb417fSDmitry Karpeev PetscMPIInt color; 15612bbb417fSDmitry Karpeev MPI_Comm comm, subcomm; 15626a4f0f73SDmitry Karpeev IS **xis = 0, **is = ois, **is_local = iis; 1563d34fcf5fSBarry Smith 1564b1a0a8a3SJed Brown PetscFunctionBegin; 15652bbb417fSDmitry Karpeev ierr = PetscObjectGetComm((PetscObject)pc, &comm);CHKERRQ(ierr); 15662bbb417fSDmitry Karpeev ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 15672bbb417fSDmitry Karpeev ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 15682bbb417fSDmitry Karpeev ierr = MatGetOwnershipRange(pc->pmat, &first, &last);CHKERRQ(ierr); 1569d34fcf5fSBarry 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) " 15702bbb417fSDmitry Karpeev "does not respect the number of degrees of freedom per grid point %D", first, last, dof); 1571d34fcf5fSBarry Smith 1572af538404SDmitry Karpeev /* Determine the number of domains with nonzero intersections with the local ownership range. */ 15732bbb417fSDmitry Karpeev s = 0; 1574af538404SDmitry Karpeev ystart = 0; 1575af538404SDmitry Karpeev for (j=0; j<Ndomains; ++j) { 1576af538404SDmitry Karpeev maxheight = N/Ndomains + ((N % Ndomains) > j); /* Maximal height of subdomain */ 1577af538404SDmitry Karpeev if (maxheight < 2) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many %D subdomains in the vertical directon for mesh height %D", Ndomains, N); 1578eec7e3faSDmitry Karpeev /* Vertical domain limits with an overlap. */ 1579eec7e3faSDmitry Karpeev ylow = PetscMax(ystart - overlap,0); 1580eec7e3faSDmitry Karpeev yhigh = PetscMin(ystart + maxheight + overlap,N); 1581b1a0a8a3SJed Brown xstart = 0; 1582af538404SDmitry Karpeev for (i=0; i<Mdomains; ++i) { 1583af538404SDmitry Karpeev maxwidth = M/Mdomains + ((M % Mdomains) > i); /* Maximal width of subdomain */ 1584af538404SDmitry 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); 1585eec7e3faSDmitry Karpeev /* Horizontal domain limits with an overlap. */ 1586eec7e3faSDmitry Karpeev xleft = PetscMax(xstart - overlap,0); 1587eec7e3faSDmitry Karpeev xright = PetscMin(xstart + maxwidth + overlap,M); 1588af538404SDmitry Karpeev /* 1589af538404SDmitry Karpeev Determine whether this subdomain intersects this processor's ownership range of pc->pmat. 1590af538404SDmitry Karpeev */ 1591c3518bceSDmitry Karpeev PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,(&xleft_loc),(&ylow_loc),(&xright_loc),(&yhigh_loc),(&nidx)); 15922fa5cd67SKarl Rupp if (nidx) ++s; 1593af538404SDmitry Karpeev xstart += maxwidth; 1594af538404SDmitry Karpeev } /* for (i = 0; i < Mdomains; ++i) */ 1595af538404SDmitry Karpeev ystart += maxheight; 1596af538404SDmitry Karpeev } /* for (j = 0; j < Ndomains; ++j) */ 15972fa5cd67SKarl Rupp 1598af538404SDmitry Karpeev /* Now we can allocate the necessary number of ISs. */ 1599af538404SDmitry Karpeev *nsub = s; 1600854ce69bSBarry Smith ierr = PetscMalloc1(*nsub,is);CHKERRQ(ierr); 1601854ce69bSBarry Smith ierr = PetscMalloc1(*nsub,is_local);CHKERRQ(ierr); 1602af538404SDmitry Karpeev s = 0; 1603af538404SDmitry Karpeev ystart = 0; 1604af538404SDmitry Karpeev for (j=0; j<Ndomains; ++j) { 1605af538404SDmitry Karpeev maxheight = N/Ndomains + ((N % Ndomains) > j); /* Maximal height of subdomain */ 1606af538404SDmitry 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); 1607af538404SDmitry Karpeev /* Vertical domain limits with an overlap. */ 1608eec7e3faSDmitry Karpeev y[0][0] = PetscMax(ystart - overlap,0); 1609eec7e3faSDmitry Karpeev y[0][1] = PetscMin(ystart + maxheight + overlap,N); 1610eec7e3faSDmitry Karpeev /* Vertical domain limits without an overlap. */ 1611eec7e3faSDmitry Karpeev y[1][0] = ystart; 1612eec7e3faSDmitry Karpeev y[1][1] = PetscMin(ystart + maxheight,N); 1613eec7e3faSDmitry Karpeev xstart = 0; 1614af538404SDmitry Karpeev for (i=0; i<Mdomains; ++i) { 1615af538404SDmitry Karpeev maxwidth = M/Mdomains + ((M % Mdomains) > i); /* Maximal width of subdomain */ 1616af538404SDmitry 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); 1617af538404SDmitry Karpeev /* Horizontal domain limits with an overlap. */ 1618eec7e3faSDmitry Karpeev x[0][0] = PetscMax(xstart - overlap,0); 1619eec7e3faSDmitry Karpeev x[0][1] = PetscMin(xstart + maxwidth + overlap,M); 1620eec7e3faSDmitry Karpeev /* Horizontal domain limits without an overlap. */ 1621eec7e3faSDmitry Karpeev x[1][0] = xstart; 1622eec7e3faSDmitry Karpeev x[1][1] = PetscMin(xstart+maxwidth,M); 16232bbb417fSDmitry Karpeev /* 16242bbb417fSDmitry Karpeev Determine whether this domain intersects this processor's ownership range of pc->pmat. 16252bbb417fSDmitry Karpeev Do this twice: first for the domains with overlaps, and once without. 16262bbb417fSDmitry Karpeev During the first pass create the subcommunicators, and use them on the second pass as well. 16272bbb417fSDmitry Karpeev */ 16282bbb417fSDmitry Karpeev for (q = 0; q < 2; ++q) { 16292bbb417fSDmitry Karpeev /* 16302bbb417fSDmitry Karpeev domain limits, (xleft, xright) and (ylow, yheigh) are adjusted 16312bbb417fSDmitry Karpeev according to whether the domain with an overlap or without is considered. 16322bbb417fSDmitry Karpeev */ 16332bbb417fSDmitry Karpeev xleft = x[q][0]; xright = x[q][1]; 16342bbb417fSDmitry Karpeev ylow = y[q][0]; yhigh = y[q][1]; 1635c3518bceSDmitry Karpeev PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,(&xleft_loc),(&ylow_loc),(&xright_loc),(&yhigh_loc),(&nidx)); 1636af538404SDmitry Karpeev nidx *= dof; 1637eec7e3faSDmitry Karpeev n[q] = nidx; 16382bbb417fSDmitry Karpeev /* 1639eec7e3faSDmitry Karpeev Based on the counted number of indices in the local domain *with an overlap*, 1640af538404SDmitry Karpeev construct a subcommunicator of all the processors supporting this domain. 1641eec7e3faSDmitry Karpeev Observe that a domain with an overlap might have nontrivial local support, 1642eec7e3faSDmitry Karpeev while the domain without an overlap might not. Hence, the decision to participate 1643eec7e3faSDmitry Karpeev in the subcommunicator must be based on the domain with an overlap. 16442bbb417fSDmitry Karpeev */ 16452bbb417fSDmitry Karpeev if (q == 0) { 16462fa5cd67SKarl Rupp if (nidx) color = 1; 16472fa5cd67SKarl Rupp else color = MPI_UNDEFINED; 16482fa5cd67SKarl Rupp 16492bbb417fSDmitry Karpeev ierr = MPI_Comm_split(comm, color, rank, &subcomm);CHKERRQ(ierr); 16502bbb417fSDmitry Karpeev } 1651af538404SDmitry Karpeev /* 1652eec7e3faSDmitry Karpeev Proceed only if the number of local indices *with an overlap* is nonzero. 1653af538404SDmitry Karpeev */ 1654eec7e3faSDmitry Karpeev if (n[0]) { 16552fa5cd67SKarl Rupp if (q == 0) xis = is; 1656af538404SDmitry Karpeev if (q == 1) { 1657af538404SDmitry Karpeev /* 1658eec7e3faSDmitry Karpeev The IS for the no-overlap subdomain shares a communicator with the overlapping domain. 1659af538404SDmitry Karpeev Moreover, if the overlap is zero, the two ISs are identical. 1660af538404SDmitry Karpeev */ 1661b1a0a8a3SJed Brown if (overlap == 0) { 1662eec7e3faSDmitry Karpeev (*is_local)[s] = (*is)[s]; 16632bbb417fSDmitry Karpeev ierr = PetscObjectReference((PetscObject)(*is)[s]);CHKERRQ(ierr); 16642bbb417fSDmitry Karpeev continue; 1665d34fcf5fSBarry Smith } else { 16666a4f0f73SDmitry Karpeev xis = is_local; 1667eec7e3faSDmitry Karpeev subcomm = ((PetscObject)(*is)[s])->comm; 16682bbb417fSDmitry Karpeev } 1669af538404SDmitry Karpeev } /* if (q == 1) */ 16700298fd71SBarry Smith idx = NULL; 1671785e854fSJed Brown ierr = PetscMalloc1(nidx,&idx);CHKERRQ(ierr); 1672eec7e3faSDmitry Karpeev if (nidx) { 16732bbb417fSDmitry Karpeev k = 0; 16742bbb417fSDmitry Karpeev for (jj=ylow_loc; jj<yhigh_loc; ++jj) { 1675af538404SDmitry Karpeev PetscInt x0 = (jj==ylow_loc) ? xleft_loc : xleft; 1676af538404SDmitry Karpeev PetscInt x1 = (jj==yhigh_loc-1) ? xright_loc : xright; 1677af538404SDmitry Karpeev kk = dof*(M*jj + x0); 16782bbb417fSDmitry Karpeev for (ii=x0; ii<x1; ++ii) { 16792bbb417fSDmitry Karpeev for (d = 0; d < dof; ++d) { 16802bbb417fSDmitry Karpeev idx[k++] = kk++; 1681b1a0a8a3SJed Brown } 1682b1a0a8a3SJed Brown } 1683b1a0a8a3SJed Brown } 1684eec7e3faSDmitry Karpeev } 16856a4f0f73SDmitry Karpeev ierr = ISCreateGeneral(subcomm,nidx,idx,PETSC_OWN_POINTER,(*xis)+s);CHKERRQ(ierr); 1686eec7e3faSDmitry Karpeev }/* if (n[0]) */ 16872bbb417fSDmitry Karpeev }/* for (q = 0; q < 2; ++q) */ 16882fa5cd67SKarl Rupp if (n[0]) ++s; 1689af538404SDmitry Karpeev xstart += maxwidth; 1690af538404SDmitry Karpeev } /* for (i = 0; i < Mdomains; ++i) */ 16912bbb417fSDmitry Karpeev ystart += maxheight; 1692af538404SDmitry Karpeev } /* for (j = 0; j < Ndomains; ++j) */ 1693b1a0a8a3SJed Brown PetscFunctionReturn(0); 1694b1a0a8a3SJed Brown } 1695b1a0a8a3SJed Brown 1696b1a0a8a3SJed Brown #undef __FUNCT__ 169706b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMGetSubdomains" 1698b1a0a8a3SJed Brown /*@C 169906b43e7eSDmitry Karpeev PCGASMGetSubdomains - Gets the subdomains supported on this processor 170006b43e7eSDmitry Karpeev for the additive Schwarz preconditioner. 1701b1a0a8a3SJed Brown 1702b1a0a8a3SJed Brown Not Collective 1703b1a0a8a3SJed Brown 1704b1a0a8a3SJed Brown Input Parameter: 1705b1a0a8a3SJed Brown . pc - the preconditioner context 1706b1a0a8a3SJed Brown 1707b1a0a8a3SJed Brown Output Parameters: 1708b1a0a8a3SJed Brown + n - the number of subdomains for this processor (default value = 1) 17090298fd71SBarry Smith . iis - the index sets that define the inner subdomains (without overlap) supported on this processor (can be NULL) 17100298fd71SBarry Smith - ois - the index sets that define the outer subdomains (with overlap) supported on this processor (can be NULL) 1711b1a0a8a3SJed Brown 1712b1a0a8a3SJed Brown 1713b1a0a8a3SJed Brown Notes: 17146a4f0f73SDmitry Karpeev The user is responsible for destroying the ISs and freeing the returned arrays. 1715b1a0a8a3SJed Brown The IS numbering is in the parallel, global numbering of the vector. 1716b1a0a8a3SJed Brown 1717b1a0a8a3SJed Brown Level: advanced 1718b1a0a8a3SJed Brown 171906b43e7eSDmitry Karpeev .keywords: PC, GASM, get, subdomains, additive Schwarz 1720b1a0a8a3SJed Brown 17218f3b4b4dSDmitry Karpeev .seealso: PCGASMSetOverlap(), PCGASMGetSubKSP(), PCGASMCreateSubdomains2D(), 17228f3b4b4dSDmitry Karpeev PCGASMSetSubdomains(), PCGASMGetSubmatrices() 1723b1a0a8a3SJed Brown @*/ 17246a4f0f73SDmitry Karpeev PetscErrorCode PCGASMGetSubdomains(PC pc,PetscInt *n,IS *iis[],IS *ois[]) 1725b1a0a8a3SJed Brown { 1726f746d493SDmitry Karpeev PC_GASM *osm; 1727b1a0a8a3SJed Brown PetscErrorCode ierr; 1728b1a0a8a3SJed Brown PetscBool match; 17296a4f0f73SDmitry Karpeev PetscInt i; 17305fd66863SKarl Rupp 1731b1a0a8a3SJed Brown PetscFunctionBegin; 1732b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1733251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match);CHKERRQ(ierr); 1734ce94432eSBarry Smith if (!match) SETERRQ2(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Incorrect object type: expected %s, got %s instead", PCGASM, ((PetscObject)pc)->type_name); 1735f746d493SDmitry Karpeev osm = (PC_GASM*)pc->data; 1736ab3e923fSDmitry Karpeev if (n) *n = osm->n; 17376a4f0f73SDmitry Karpeev if (iis) { 1738785e854fSJed Brown ierr = PetscMalloc1(osm->n, iis);CHKERRQ(ierr); 17396a4f0f73SDmitry Karpeev } 17406a4f0f73SDmitry Karpeev if (ois) { 1741785e854fSJed Brown ierr = PetscMalloc1(osm->n, ois);CHKERRQ(ierr); 17426a4f0f73SDmitry Karpeev } 17436a4f0f73SDmitry Karpeev if (iis || ois) { 17446a4f0f73SDmitry Karpeev for (i = 0; i < osm->n; ++i) { 17456a4f0f73SDmitry Karpeev if (iis) (*iis)[i] = osm->iis[i]; 17466a4f0f73SDmitry Karpeev if (ois) (*ois)[i] = osm->ois[i]; 17476a4f0f73SDmitry Karpeev } 1748b1a0a8a3SJed Brown } 1749b1a0a8a3SJed Brown PetscFunctionReturn(0); 1750b1a0a8a3SJed Brown } 1751b1a0a8a3SJed Brown 1752b1a0a8a3SJed Brown #undef __FUNCT__ 175306b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMGetSubmatrices" 1754b1a0a8a3SJed Brown /*@C 175506b43e7eSDmitry Karpeev PCGASMGetSubmatrices - Gets the local submatrices (for this processor 1756b1a0a8a3SJed Brown only) for the additive Schwarz preconditioner. 1757b1a0a8a3SJed Brown 1758b1a0a8a3SJed Brown Not Collective 1759b1a0a8a3SJed Brown 1760b1a0a8a3SJed Brown Input Parameter: 1761b1a0a8a3SJed Brown . pc - the preconditioner context 1762b1a0a8a3SJed Brown 1763b1a0a8a3SJed Brown Output Parameters: 1764b1a0a8a3SJed Brown + n - the number of matrices for this processor (default value = 1) 1765b1a0a8a3SJed Brown - mat - the matrices 1766b1a0a8a3SJed Brown 176706b43e7eSDmitry Karpeev Notes: matrices returned by this routine have the same communicators as the index sets (IS) 17688f3b4b4dSDmitry Karpeev used to define subdomains in PCGASMSetSubdomains() 1769b1a0a8a3SJed Brown Level: advanced 1770b1a0a8a3SJed Brown 1771f746d493SDmitry Karpeev .keywords: PC, GASM, set, local, subdomains, additive Schwarz, block Jacobi 1772b1a0a8a3SJed Brown 17738f3b4b4dSDmitry Karpeev .seealso: PCGASMSetOverlap(), PCGASMGetSubKSP(), 177406b43e7eSDmitry Karpeev PCGASMCreateSubdomains2D(), PCGASMSetSubdomains(), PCGASMGetSubdomains() 1775b1a0a8a3SJed Brown @*/ 177606b43e7eSDmitry Karpeev PetscErrorCode PCGASMGetSubmatrices(PC pc,PetscInt *n,Mat *mat[]) 1777b1a0a8a3SJed Brown { 1778f746d493SDmitry Karpeev PC_GASM *osm; 1779b1a0a8a3SJed Brown PetscErrorCode ierr; 1780b1a0a8a3SJed Brown PetscBool match; 1781b1a0a8a3SJed Brown 1782b1a0a8a3SJed Brown PetscFunctionBegin; 1783b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1784b1a0a8a3SJed Brown PetscValidIntPointer(n,2); 1785b1a0a8a3SJed Brown if (mat) PetscValidPointer(mat,3); 1786ce94432eSBarry Smith if (!pc->setupcalled) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call after KSPSetUP() or PCSetUp()."); 1787251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match);CHKERRQ(ierr); 1788ce94432eSBarry Smith if (!match) SETERRQ2(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Expected %s, got %s instead", PCGASM, ((PetscObject)pc)->type_name); 1789f746d493SDmitry Karpeev osm = (PC_GASM*)pc->data; 1790ab3e923fSDmitry Karpeev if (n) *n = osm->n; 1791b1a0a8a3SJed Brown if (mat) *mat = osm->pmat; 1792b1a0a8a3SJed Brown PetscFunctionReturn(0); 1793b1a0a8a3SJed Brown } 1794d709ea83SDmitry Karpeev 1795d709ea83SDmitry Karpeev #undef __FUNCT__ 17968f3b4b4dSDmitry Karpeev #define __FUNCT__ "PCGASMSetUseDMSubdomains" 1797d709ea83SDmitry Karpeev /*@ 17988f3b4b4dSDmitry Karpeev PCGASMSetUseDMSubdomains - Indicates whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible. 1799d709ea83SDmitry Karpeev Logically Collective 1800d709ea83SDmitry Karpeev 1801d709ea83SDmitry Karpeev Input Parameter: 1802d709ea83SDmitry Karpeev + pc - the preconditioner 1803d709ea83SDmitry Karpeev - flg - boolean indicating whether to use subdomains defined by the DM 1804d709ea83SDmitry Karpeev 1805d709ea83SDmitry Karpeev Options Database Key: 18068f3b4b4dSDmitry Karpeev . -pc_gasm_dm_subdomains -pc_gasm_overlap -pc_gasm_total_subdomains 1807d709ea83SDmitry Karpeev 1808d709ea83SDmitry Karpeev Level: intermediate 1809d709ea83SDmitry Karpeev 1810d709ea83SDmitry Karpeev Notes: 18118f3b4b4dSDmitry Karpeev PCGASMSetSubdomains(), PCGASMSetTotalSubdomains() or PCGASMSetOverlap() take precedence over PCGASMSetUseDMSubdomains(), 18128f3b4b4dSDmitry Karpeev so setting PCGASMSetSubdomains() with nontrivial subdomain ISs or any of PCGASMSetTotalSubdomains() and PCGASMSetOverlap() 18138f3b4b4dSDmitry Karpeev automatically turns the latter off. 1814d709ea83SDmitry Karpeev 1815d709ea83SDmitry Karpeev .keywords: PC, ASM, DM, set, subdomains, additive Schwarz 1816d709ea83SDmitry Karpeev 18178f3b4b4dSDmitry Karpeev .seealso: PCGASMGetUseDMSubdomains(), PCGASMSetSubdomains(), PCGASMSetOverlap() 1818d709ea83SDmitry Karpeev PCGASMCreateSubdomains2D() 1819d709ea83SDmitry Karpeev @*/ 18208f3b4b4dSDmitry Karpeev PetscErrorCode PCGASMSetUseDMSubdomains(PC pc,PetscBool flg) 1821d709ea83SDmitry Karpeev { 1822d709ea83SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 1823d709ea83SDmitry Karpeev PetscErrorCode ierr; 1824d709ea83SDmitry Karpeev PetscBool match; 1825d709ea83SDmitry Karpeev 1826d709ea83SDmitry Karpeev PetscFunctionBegin; 1827d709ea83SDmitry Karpeev PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1828d709ea83SDmitry Karpeev PetscValidLogicalCollectiveBool(pc,flg,2); 1829d709ea83SDmitry Karpeev if (pc->setupcalled) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Not for a setup PC."); 1830d709ea83SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match);CHKERRQ(ierr); 1831d709ea83SDmitry Karpeev if (match) { 18328f3b4b4dSDmitry Karpeev if (!osm->user_subdomains && osm->N == PETSC_DETERMINE && osm->overlap < 0) { 1833d709ea83SDmitry Karpeev osm->dm_subdomains = flg; 1834d709ea83SDmitry Karpeev } 18358f3b4b4dSDmitry Karpeev } 1836d709ea83SDmitry Karpeev PetscFunctionReturn(0); 1837d709ea83SDmitry Karpeev } 1838d709ea83SDmitry Karpeev 1839d709ea83SDmitry Karpeev #undef __FUNCT__ 18408f3b4b4dSDmitry Karpeev #define __FUNCT__ "PCGASMGetUseDMSubdomains" 1841d709ea83SDmitry Karpeev /*@ 18428f3b4b4dSDmitry Karpeev PCGASMGetUseDMSubdomains - Returns flag indicating whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible. 1843d709ea83SDmitry Karpeev Not Collective 1844d709ea83SDmitry Karpeev 1845d709ea83SDmitry Karpeev Input Parameter: 1846d709ea83SDmitry Karpeev . pc - the preconditioner 1847d709ea83SDmitry Karpeev 1848d709ea83SDmitry Karpeev Output Parameter: 1849d709ea83SDmitry Karpeev . flg - boolean indicating whether to use subdomains defined by the DM 1850d709ea83SDmitry Karpeev 1851d709ea83SDmitry Karpeev Level: intermediate 1852d709ea83SDmitry Karpeev 1853d709ea83SDmitry Karpeev .keywords: PC, ASM, DM, set, subdomains, additive Schwarz 1854d709ea83SDmitry Karpeev 18558f3b4b4dSDmitry Karpeev .seealso: PCGASMSetUseDMSubdomains(), PCGASMSetOverlap() 1856d709ea83SDmitry Karpeev PCGASMCreateSubdomains2D() 1857d709ea83SDmitry Karpeev @*/ 18588f3b4b4dSDmitry Karpeev PetscErrorCode PCGASMGetUseDMSubdomains(PC pc,PetscBool* flg) 1859d709ea83SDmitry Karpeev { 1860d709ea83SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 1861d709ea83SDmitry Karpeev PetscErrorCode ierr; 1862d709ea83SDmitry Karpeev PetscBool match; 1863d709ea83SDmitry Karpeev 1864d709ea83SDmitry Karpeev PetscFunctionBegin; 1865d709ea83SDmitry Karpeev PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1866d709ea83SDmitry Karpeev PetscValidPointer(flg,2); 1867d709ea83SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match);CHKERRQ(ierr); 1868d709ea83SDmitry Karpeev if (match) { 1869d709ea83SDmitry Karpeev if (flg) *flg = osm->dm_subdomains; 1870d709ea83SDmitry Karpeev } 1871d709ea83SDmitry Karpeev PetscFunctionReturn(0); 1872d709ea83SDmitry Karpeev } 1873