1b1a0a8a3SJed Brown /* 2f746d493SDmitry Karpeev This file defines an "generalized" additive Schwarz preconditioner for any Mat implementation. 36a4f0f73SDmitry Karpeev In this version each processor may intersect multiple subdomains and any subdomain may 46a4f0f73SDmitry Karpeev intersect multiple processors. Intersections of subdomains with processors are called *local 56a4f0f73SDmitry Karpeev subdomains*. 6b1a0a8a3SJed Brown 78f3b4b4dSDmitry Karpeev N - total number of distinct global subdomains (set explicitly in PCGASMSetTotalSubdomains() or implicitly PCGASMSetSubdomains() and then calculated in PCSetUp_GASM()) 86a4f0f73SDmitry Karpeev n - actual number of local subdomains on this processor (set in PCGASMSetSubdomains() or calculated in PCGASMSetTotalSubdomains()) 98f3b4b4dSDmitry Karpeev nmax - maximum number of local subdomains per processor (calculated in PCSetUp_GASM()) 10b1a0a8a3SJed Brown */ 11af0996ceSBarry Smith #include <petsc/private/pcimpl.h> /*I "petscpc.h" I*/ 121e25c274SJed Brown #include <petscdm.h> 13b1a0a8a3SJed Brown 14b1a0a8a3SJed Brown typedef struct { 15f746d493SDmitry Karpeev PetscInt N,n,nmax; 16b1a0a8a3SJed Brown PetscInt overlap; /* overlap requested by user */ 178f3b4b4dSDmitry Karpeev PCGASMType type; /* use reduced interpolation, restriction or both */ 188f3b4b4dSDmitry Karpeev PetscBool type_set; /* if user set this value (so won't change it for symmetric problems) */ 198f3b4b4dSDmitry Karpeev PetscBool same_subdomain_solvers; /* flag indicating whether all local solvers are same */ 208f3b4b4dSDmitry Karpeev PetscBool sort_indices; /* flag to sort subdomain indices */ 218f3b4b4dSDmitry Karpeev PetscBool user_subdomains; /* whether the user set explicit subdomain index sets -- keep them on PCReset() */ 228f3b4b4dSDmitry Karpeev PetscBool dm_subdomains; /* whether DM is allowed to define subdomains */ 23ea91fabdSFande Kong PetscBool hierarchicalpartitioning; 248f3b4b4dSDmitry Karpeev IS *ois; /* index sets that define the outer (conceptually, overlapping) subdomains */ 258f3b4b4dSDmitry Karpeev IS *iis; /* index sets that define the inner (conceptually, nonoverlapping) subdomains */ 268f3b4b4dSDmitry Karpeev KSP *ksp; /* linear solvers for each subdomain */ 278f3b4b4dSDmitry Karpeev Mat *pmat; /* subdomain block matrices */ 28f746d493SDmitry Karpeev Vec gx,gy; /* Merged work vectors */ 29f746d493SDmitry Karpeev Vec *x,*y; /* Split work vectors; storage aliases pieces of storage of the above merged vectors. */ 306a4f0f73SDmitry Karpeev VecScatter gorestriction; /* merged restriction to disjoint union of outer subdomains */ 316a4f0f73SDmitry Karpeev VecScatter girestriction; /* merged restriction to disjoint union of inner subdomains */ 32ea91fabdSFande Kong VecScatter pctoouter; 33ea91fabdSFande Kong IS permutationIS; 34ea91fabdSFande Kong Mat permutationP; 35ea91fabdSFande Kong Mat pcmat; 36ea91fabdSFande Kong Vec pcx,pcy; 37f746d493SDmitry Karpeev } PC_GASM; 38b1a0a8a3SJed Brown 398f3b4b4dSDmitry Karpeev static PetscErrorCode PCGASMComputeGlobalSubdomainNumbering_Private(PC pc,PetscInt **numbering,PetscInt **permutation) 408f3b4b4dSDmitry Karpeev { 418f3b4b4dSDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 428f3b4b4dSDmitry Karpeev PetscInt i; 438f3b4b4dSDmitry Karpeev PetscErrorCode ierr; 448f3b4b4dSDmitry Karpeev 458f3b4b4dSDmitry Karpeev PetscFunctionBegin; 468f3b4b4dSDmitry Karpeev /* Determine the number of globally-distinct subdomains and compute a global numbering for them. */ 478f3b4b4dSDmitry Karpeev ierr = PetscMalloc2(osm->n,numbering,osm->n,permutation);CHKERRQ(ierr); 488f3b4b4dSDmitry Karpeev ierr = PetscObjectsListGetGlobalNumbering(PetscObjectComm((PetscObject)pc),osm->n,(PetscObject*)osm->iis,NULL,*numbering);CHKERRQ(ierr); 498f3b4b4dSDmitry Karpeev for (i = 0; i < osm->n; ++i) (*permutation)[i] = i; 508f3b4b4dSDmitry Karpeev ierr = PetscSortIntWithPermutation(osm->n,*numbering,*permutation);CHKERRQ(ierr); 518f3b4b4dSDmitry Karpeev PetscFunctionReturn(0); 528f3b4b4dSDmitry Karpeev } 538f3b4b4dSDmitry Karpeev 5406b43e7eSDmitry Karpeev static PetscErrorCode PCGASMSubdomainView_Private(PC pc, PetscInt i, PetscViewer viewer) 55af538404SDmitry Karpeev { 56af538404SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 5706b43e7eSDmitry Karpeev PetscInt j,nidx; 58af538404SDmitry Karpeev const PetscInt *idx; 5906b43e7eSDmitry Karpeev PetscViewer sviewer; 60af538404SDmitry Karpeev char *cidx; 61af538404SDmitry Karpeev PetscErrorCode ierr; 62af538404SDmitry Karpeev 63af538404SDmitry Karpeev PetscFunctionBegin; 64ce94432eSBarry Smith if (i < 0 || i > osm->n) SETERRQ2(PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONG, "Invalid subdomain %D: must nonnegative and less than %D", i, osm->n); 654bde246dSDmitry Karpeev /* Inner subdomains. */ 664bde246dSDmitry Karpeev ierr = ISGetLocalSize(osm->iis[i], &nidx);CHKERRQ(ierr); 674bde246dSDmitry Karpeev /* 684bde246dSDmitry Karpeev No more than 15 characters per index plus a space. 694bde246dSDmitry Karpeev PetscViewerStringSPrintf requires a string of size at least 2, so use (nidx+1) instead of nidx, 704bde246dSDmitry Karpeev in case nidx == 0. That will take care of the space for the trailing '\0' as well. 714bde246dSDmitry Karpeev For nidx == 0, the whole string 16 '\0'. 724bde246dSDmitry Karpeev */ 7389d949e2SBarry Smith ierr = PetscMalloc1(16*(nidx+1)+1, &cidx);CHKERRQ(ierr); 744bde246dSDmitry Karpeev ierr = PetscViewerStringOpen(PETSC_COMM_SELF, cidx, 16*(nidx+1)+1, &sviewer);CHKERRQ(ierr); 754bde246dSDmitry Karpeev ierr = ISGetIndices(osm->iis[i], &idx);CHKERRQ(ierr); 764bde246dSDmitry Karpeev for (j = 0; j < nidx; ++j) { 774bde246dSDmitry Karpeev ierr = PetscViewerStringSPrintf(sviewer, "%D ", idx[j]);CHKERRQ(ierr); 784bde246dSDmitry Karpeev } 794bde246dSDmitry Karpeev ierr = ISRestoreIndices(osm->iis[i],&idx);CHKERRQ(ierr); 804bde246dSDmitry Karpeev ierr = PetscViewerDestroy(&sviewer);CHKERRQ(ierr); 814bde246dSDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer, "Inner subdomain:\n");CHKERRQ(ierr); 824bde246dSDmitry Karpeev ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 831575c14dSBarry Smith ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 844bde246dSDmitry Karpeev ierr = PetscViewerASCIISynchronizedPrintf(viewer, "%s", cidx);CHKERRQ(ierr); 854bde246dSDmitry Karpeev ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 861575c14dSBarry Smith ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 874bde246dSDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer, "\n");CHKERRQ(ierr); 884bde246dSDmitry Karpeev ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 894bde246dSDmitry Karpeev ierr = PetscFree(cidx);CHKERRQ(ierr); 904bde246dSDmitry Karpeev /* Outer subdomains. */ 916a4f0f73SDmitry Karpeev ierr = ISGetLocalSize(osm->ois[i], &nidx);CHKERRQ(ierr); 92eec7e3faSDmitry Karpeev /* 93eec7e3faSDmitry Karpeev No more than 15 characters per index plus a space. 94eec7e3faSDmitry Karpeev PetscViewerStringSPrintf requires a string of size at least 2, so use (nidx+1) instead of nidx, 95eec7e3faSDmitry Karpeev in case nidx == 0. That will take care of the space for the trailing '\0' as well. 96eec7e3faSDmitry Karpeev For nidx == 0, the whole string 16 '\0'. 97eec7e3faSDmitry Karpeev */ 9889d949e2SBarry Smith ierr = PetscMalloc1(16*(nidx+1)+1, &cidx);CHKERRQ(ierr); 9906b43e7eSDmitry Karpeev ierr = PetscViewerStringOpen(PETSC_COMM_SELF, cidx, 16*(nidx+1)+1, &sviewer);CHKERRQ(ierr); 1006a4f0f73SDmitry Karpeev ierr = ISGetIndices(osm->ois[i], &idx);CHKERRQ(ierr); 101af538404SDmitry Karpeev for (j = 0; j < nidx; ++j) { 102af538404SDmitry Karpeev ierr = PetscViewerStringSPrintf(sviewer,"%D ", idx[j]);CHKERRQ(ierr); 103af538404SDmitry Karpeev } 1046bf464f9SBarry Smith ierr = PetscViewerDestroy(&sviewer);CHKERRQ(ierr); 1056a4f0f73SDmitry Karpeev ierr = ISRestoreIndices(osm->ois[i],&idx);CHKERRQ(ierr); 1066a4f0f73SDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer, "Outer subdomain:\n");CHKERRQ(ierr); 10706b43e7eSDmitry Karpeev ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 1081575c14dSBarry Smith ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 109af538404SDmitry Karpeev ierr = PetscViewerASCIISynchronizedPrintf(viewer, "%s", cidx);CHKERRQ(ierr); 110af538404SDmitry Karpeev ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 1111575c14dSBarry Smith ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 112af538404SDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer, "\n");CHKERRQ(ierr); 11306b43e7eSDmitry Karpeev ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 114af538404SDmitry Karpeev ierr = PetscFree(cidx);CHKERRQ(ierr); 11506b43e7eSDmitry Karpeev PetscFunctionReturn(0); 11606b43e7eSDmitry Karpeev } 11706b43e7eSDmitry Karpeev 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; 132c5929fdfSBarry Smith ierr = PetscOptionsGetBool(NULL,prefix,"-pc_gasm_print_subdomains",&doprint,NULL);CHKERRQ(ierr); 1338f3b4b4dSDmitry Karpeev if (!doprint) PetscFunctionReturn(0); 134c5929fdfSBarry Smith ierr = PetscOptionsGetString(NULL,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) { 1493f08860eSBarry Smith ierr = PetscViewerGetSubViewer(viewer,((PetscObject)osm->ois[d])->comm, &sviewer);CHKERRQ(ierr); 1504bde246dSDmitry Karpeev ierr = PCGASMSubdomainView_Private(pc,d,sviewer);CHKERRQ(ierr); 1513f08860eSBarry Smith 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 163f746d493SDmitry Karpeev static PetscErrorCode PCView_GASM(PC pc,PetscViewer viewer) 164b1a0a8a3SJed Brown { 165f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 166af538404SDmitry Karpeev const char *prefix; 167b1a0a8a3SJed Brown PetscErrorCode ierr; 168af538404SDmitry Karpeev PetscMPIInt rank, size; 1698f3b4b4dSDmitry Karpeev PetscInt bsz; 17006b43e7eSDmitry Karpeev PetscBool iascii,view_subdomains=PETSC_FALSE; 171b1a0a8a3SJed Brown PetscViewer sviewer; 1728f3b4b4dSDmitry Karpeev PetscInt count, l; 1736a4f0f73SDmitry Karpeev char overlap[256] = "user-defined overlap"; 1746a4f0f73SDmitry Karpeev char gsubdomains[256] = "unknown total number of subdomains"; 17506b43e7eSDmitry Karpeev char msubdomains[256] = "unknown max number of local subdomains"; 1768f3b4b4dSDmitry Karpeev PetscInt *numbering,*permutation;/* global numbering of locally-supported subdomains and the permutation from the local ordering */ 17711aeaf0aSBarry Smith 178b1a0a8a3SJed Brown PetscFunctionBegin; 179ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size);CHKERRQ(ierr); 180ce94432eSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank);CHKERRQ(ierr); 18106b43e7eSDmitry Karpeev 18206b43e7eSDmitry Karpeev if (osm->overlap >= 0) { 18306b43e7eSDmitry Karpeev ierr = PetscSNPrintf(overlap,sizeof(overlap),"requested amount of overlap = %D",osm->overlap);CHKERRQ(ierr); 18406b43e7eSDmitry Karpeev } 1858f3b4b4dSDmitry Karpeev if (osm->N != PETSC_DETERMINE) { 1868f3b4b4dSDmitry Karpeev ierr = PetscSNPrintf(gsubdomains, sizeof(gsubdomains), "total number of subdomains = %D",osm->N);CHKERRQ(ierr); 18706b43e7eSDmitry Karpeev } 1888f3b4b4dSDmitry Karpeev if (osm->nmax != PETSC_DETERMINE) { 18906b43e7eSDmitry Karpeev ierr = PetscSNPrintf(msubdomains,sizeof(msubdomains),"max number of local subdomains = %D",osm->nmax);CHKERRQ(ierr); 19006b43e7eSDmitry Karpeev } 19106b43e7eSDmitry Karpeev 192af538404SDmitry Karpeev ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 193c5929fdfSBarry Smith ierr = PetscOptionsGetBool(NULL,prefix,"-pc_gasm_view_subdomains",&view_subdomains,NULL);CHKERRQ(ierr); 19406b43e7eSDmitry Karpeev 19506b43e7eSDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 196b1a0a8a3SJed Brown if (iascii) { 19706b43e7eSDmitry Karpeev /* 19806b43e7eSDmitry Karpeev Make sure the viewer has a name. Otherwise this may cause a deadlock when creating a subcomm viewer: 19906b43e7eSDmitry Karpeev the subcomm viewer will attempt to inherit the viewer's name, which, if not set, will be constructed 20006b43e7eSDmitry Karpeev collectively on the comm. 20106b43e7eSDmitry Karpeev */ 20206b43e7eSDmitry Karpeev ierr = PetscObjectName((PetscObject)viewer);CHKERRQ(ierr); 20306b43e7eSDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer," Restriction/interpolation type: %s\n",PCGASMTypes[osm->type]);CHKERRQ(ierr); 20406b43e7eSDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer," %s\n",overlap);CHKERRQ(ierr); 20506b43e7eSDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer," %s\n",gsubdomains);CHKERRQ(ierr); 20606b43e7eSDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer," %s\n",msubdomains);CHKERRQ(ierr); 2071575c14dSBarry Smith ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 2088f3b4b4dSDmitry Karpeev ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d|%d] number of locally-supported subdomains = %D\n",rank,size,osm->n);CHKERRQ(ierr); 209af538404SDmitry Karpeev ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 2101575c14dSBarry Smith ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 2116a4f0f73SDmitry Karpeev /* Cannot take advantage of osm->same_subdomain_solvers without a global numbering of subdomains. */ 21206b43e7eSDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer," Subdomain solver info is as follows:\n");CHKERRQ(ierr); 213b1a0a8a3SJed Brown ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 214b1a0a8a3SJed Brown ierr = PetscViewerASCIIPrintf(viewer," - - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr); 21506b43e7eSDmitry Karpeev /* Make sure that everybody waits for the banner to be printed. */ 216ce94432eSBarry Smith ierr = MPI_Barrier(PetscObjectComm((PetscObject)viewer));CHKERRQ(ierr); 2174bde246dSDmitry Karpeev /* Now let subdomains go one at a time in the global numbering order and print their subdomain/solver info. */ 2188f3b4b4dSDmitry Karpeev ierr = PCGASMComputeGlobalSubdomainNumbering_Private(pc,&numbering,&permutation);CHKERRQ(ierr); 21906b43e7eSDmitry Karpeev l = 0; 2208f3b4b4dSDmitry Karpeev for (count = 0; count < osm->N; ++count) { 22106b43e7eSDmitry Karpeev PetscMPIInt srank, ssize; 22206b43e7eSDmitry Karpeev if (l<osm->n) { 22306b43e7eSDmitry Karpeev PetscInt d = permutation[l]; /* d is the local number of the l-th smallest (in the global ordering) among the locally supported subdomains */ 22406b43e7eSDmitry Karpeev if (numbering[d] == count) { 2256a4f0f73SDmitry Karpeev ierr = MPI_Comm_size(((PetscObject)osm->ois[d])->comm, &ssize);CHKERRQ(ierr); 2266a4f0f73SDmitry Karpeev ierr = MPI_Comm_rank(((PetscObject)osm->ois[d])->comm, &srank);CHKERRQ(ierr); 2273f08860eSBarry Smith ierr = PetscViewerGetSubViewer(viewer,((PetscObject)osm->ois[d])->comm, &sviewer);CHKERRQ(ierr); 2286a4f0f73SDmitry Karpeev ierr = ISGetLocalSize(osm->ois[d],&bsz);CHKERRQ(ierr); 2298f3b4b4dSDmitry 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); 2306a4f0f73SDmitry Karpeev ierr = PetscViewerFlush(sviewer);CHKERRQ(ierr); 23106b43e7eSDmitry Karpeev if (view_subdomains) { 23206b43e7eSDmitry Karpeev ierr = PCGASMSubdomainView_Private(pc,d,sviewer);CHKERRQ(ierr); 23306b43e7eSDmitry Karpeev } 2346a4f0f73SDmitry Karpeev if (!pc->setupcalled) { 2351575c14dSBarry Smith ierr = PetscViewerASCIIPrintf(sviewer, " Solver not set up yet: PCSetUp() not yet called\n");CHKERRQ(ierr); 2362fa5cd67SKarl Rupp } else { 23706b43e7eSDmitry Karpeev ierr = KSPView(osm->ksp[d],sviewer);CHKERRQ(ierr); 2386a4f0f73SDmitry Karpeev } 23906b43e7eSDmitry Karpeev ierr = PetscViewerASCIIPrintf(sviewer," - - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr); 2407b23a99aSBarry Smith ierr = PetscViewerFlush(sviewer);CHKERRQ(ierr); 2413f08860eSBarry Smith ierr = PetscViewerRestoreSubViewer(viewer,((PetscObject)osm->ois[d])->comm, &sviewer);CHKERRQ(ierr); 24206b43e7eSDmitry Karpeev ++l; 243b1a0a8a3SJed Brown } 24406b43e7eSDmitry Karpeev } 245ce94432eSBarry Smith ierr = MPI_Barrier(PetscObjectComm((PetscObject)pc));CHKERRQ(ierr); 246b1a0a8a3SJed Brown } 2478f3b4b4dSDmitry Karpeev ierr = PetscFree2(numbering,permutation);CHKERRQ(ierr); 248b1a0a8a3SJed Brown ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 249b1a0a8a3SJed Brown } 250c5e4d11fSDmitry Karpeev ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 2511575c14dSBarry Smith ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 252b1a0a8a3SJed Brown PetscFunctionReturn(0); 253b1a0a8a3SJed Brown } 254b1a0a8a3SJed Brown 2558f3b4b4dSDmitry Karpeev PETSC_INTERN PetscErrorCode PCGASMCreateLocalSubdomains(Mat A, PetscInt nloc, IS *iis[]); 256af538404SDmitry Karpeev 257ea91fabdSFande Kong 258ea91fabdSFande Kong 259ea91fabdSFande Kong PetscErrorCode PCGASMSetHierarchicalPartitioning(PC pc) 260ea91fabdSFande Kong { 261ea91fabdSFande Kong PC_GASM *osm = (PC_GASM*)pc->data; 262ea91fabdSFande Kong MatPartitioning part; 263ea91fabdSFande Kong MPI_Comm comm; 264ea91fabdSFande Kong PetscMPIInt size; 265ea91fabdSFande Kong PetscInt nlocalsubdomains,fromrows_localsize; 266ea91fabdSFande Kong IS partitioning,fromrows,isn; 267ea91fabdSFande Kong Vec outervec; 268ea91fabdSFande Kong PetscErrorCode ierr; 269ea91fabdSFande Kong 270ea91fabdSFande Kong PetscFunctionBegin; 271ea91fabdSFande Kong ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 272ea91fabdSFande Kong ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 273ea91fabdSFande Kong /* we do not need a hierarchical partitioning when 274ea91fabdSFande Kong * the total number of subdomains is consistent with 275ea91fabdSFande Kong * the number of MPI tasks. 276ea91fabdSFande Kong * For the following cases, we do not need to use HP 277ea91fabdSFande Kong * */ 278670417b2SFande Kong if(osm->N==PETSC_DETERMINE || osm->N>=size || osm->N==1) PetscFunctionReturn(0); 279670417b2SFande Kong if(size%osm->N != 0) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_ARG_INCOMP,"have to specify the total number of subdomains %D to be a factor of the number of processors %d \n",osm->N,size); 280ea91fabdSFande Kong nlocalsubdomains = size/osm->N; 281ea91fabdSFande Kong osm->n = 1; 282ea91fabdSFande Kong ierr = MatPartitioningCreate(comm,&part);CHKERRQ(ierr); 283ea91fabdSFande Kong ierr = MatPartitioningSetAdjacency(part,pc->pmat);CHKERRQ(ierr); 284ea91fabdSFande Kong ierr = MatPartitioningSetType(part,MATPARTITIONINGHIERARCH);CHKERRQ(ierr); 285ea91fabdSFande Kong ierr = MatPartitioningHierarchicalSetNcoarseparts(part,osm->N);CHKERRQ(ierr); 286ea91fabdSFande Kong ierr = MatPartitioningHierarchicalSetNfineparts(part,nlocalsubdomains);CHKERRQ(ierr); 287ea91fabdSFande Kong ierr = MatPartitioningSetFromOptions(part);CHKERRQ(ierr); 288ea91fabdSFande Kong /* get new processor owner number of each vertex */ 289ea91fabdSFande Kong ierr = MatPartitioningApply(part,&partitioning);CHKERRQ(ierr); 290ea91fabdSFande Kong ierr = ISBuildTwoSided(partitioning,NULL,&fromrows);CHKERRQ(ierr); 291ea91fabdSFande Kong ierr = ISPartitioningToNumbering(partitioning,&isn);CHKERRQ(ierr); 292ea91fabdSFande Kong ierr = ISDestroy(&isn);CHKERRQ(ierr); 293ea91fabdSFande Kong ierr = ISGetLocalSize(fromrows,&fromrows_localsize);CHKERRQ(ierr); 294ea91fabdSFande Kong ierr = MatPartitioningDestroy(&part);CHKERRQ(ierr); 295ea91fabdSFande Kong ierr = MatCreateVecs(pc->pmat,&outervec,NULL);CHKERRQ(ierr); 296ea91fabdSFande Kong ierr = VecCreateMPI(comm,fromrows_localsize,PETSC_DETERMINE,&(osm->pcx));CHKERRQ(ierr); 297ea91fabdSFande Kong ierr = VecDuplicate(osm->pcx,&(osm->pcy));CHKERRQ(ierr); 298*9448b7f1SJunchao Zhang ierr = VecScatterCreate(osm->pcx,NULL,outervec,fromrows,&(osm->pctoouter));CHKERRQ(ierr); 2997dae84e0SHong Zhang ierr = MatCreateSubMatrix(pc->pmat,fromrows,fromrows,MAT_INITIAL_MATRIX,&(osm->permutationP));CHKERRQ(ierr); 300ea91fabdSFande Kong ierr = PetscObjectReference((PetscObject)fromrows);CHKERRQ(ierr); 301ea91fabdSFande Kong osm->permutationIS = fromrows; 302ea91fabdSFande Kong osm->pcmat = pc->pmat; 303ea91fabdSFande Kong ierr = PetscObjectReference((PetscObject)osm->permutationP);CHKERRQ(ierr); 304ea91fabdSFande Kong pc->pmat = osm->permutationP; 305ea91fabdSFande Kong ierr = VecDestroy(&outervec);CHKERRQ(ierr); 306ea91fabdSFande Kong ierr = ISDestroy(&fromrows);CHKERRQ(ierr); 307ea91fabdSFande Kong ierr = ISDestroy(&partitioning);CHKERRQ(ierr); 308ea91fabdSFande Kong osm->n = PETSC_DETERMINE; 309ea91fabdSFande Kong PetscFunctionReturn(0); 310ea91fabdSFande Kong } 311ea91fabdSFande Kong 312ea91fabdSFande Kong 313ea91fabdSFande Kong 314f746d493SDmitry Karpeev static PetscErrorCode PCSetUp_GASM(PC pc) 315b1a0a8a3SJed Brown { 316f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 317b1a0a8a3SJed Brown PetscErrorCode ierr; 318b1a0a8a3SJed Brown PetscBool symset,flg; 319ea91fabdSFande Kong PetscInt i,nInnerIndices,nTotalInnerIndices; 320c10bc1a1SDmitry Karpeev PetscMPIInt rank, size; 321b1a0a8a3SJed Brown MatReuse scall = MAT_REUSE_MATRIX; 322b1a0a8a3SJed Brown KSP ksp; 323b1a0a8a3SJed Brown PC subpc; 324b1a0a8a3SJed Brown const char *prefix,*pprefix; 325f746d493SDmitry Karpeev Vec x,y; 3266a4f0f73SDmitry Karpeev PetscInt oni; /* Number of indices in the i-th local outer subdomain. */ 3276a4f0f73SDmitry Karpeev const PetscInt *oidxi; /* Indices from the i-th subdomain local outer subdomain. */ 3286a4f0f73SDmitry Karpeev PetscInt on; /* Number of indices in the disjoint union of local outer subdomains. */ 3296a4f0f73SDmitry Karpeev PetscInt *oidx; /* Indices in the disjoint union of local outer subdomains. */ 3306a4f0f73SDmitry Karpeev IS gois; /* Disjoint union the global indices of outer subdomains. */ 3316a4f0f73SDmitry Karpeev IS goid; /* Identity IS of the size of the disjoint union of outer subdomains. */ 332f746d493SDmitry Karpeev PetscScalar *gxarray, *gyarray; 333930d09c1SFande Kong PetscInt gostart; /* Start of locally-owned indices in the vectors -- osm->gx,osm->gy -- over the disjoint union of outer subdomains. */ 3348f3b4b4dSDmitry Karpeev PetscInt num_subdomains = 0; 3350298fd71SBarry Smith DM *subdomain_dm = NULL; 3368f3b4b4dSDmitry Karpeev char **subdomain_names = NULL; 337f771a274SFande Kong PetscInt *numbering; 3388f3b4b4dSDmitry Karpeev 339b1a0a8a3SJed Brown 340b1a0a8a3SJed Brown PetscFunctionBegin; 341ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr); 342ce94432eSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);CHKERRQ(ierr); 343b1a0a8a3SJed Brown if (!pc->setupcalled) { 344ea91fabdSFande Kong /* use a hierarchical partitioning */ 345ea91fabdSFande Kong if(osm->hierarchicalpartitioning){ 346ea91fabdSFande Kong ierr = PCGASMSetHierarchicalPartitioning(pc);CHKERRQ(ierr); 347ea91fabdSFande Kong } 348b1a0a8a3SJed Brown if (!osm->type_set) { 349b1a0a8a3SJed Brown ierr = MatIsSymmetricKnown(pc->pmat,&symset,&flg);CHKERRQ(ierr); 3502fa5cd67SKarl Rupp if (symset && flg) osm->type = PC_GASM_BASIC; 351b1a0a8a3SJed Brown } 352b1a0a8a3SJed Brown 3538f3b4b4dSDmitry Karpeev if (osm->n == PETSC_DETERMINE) { 3548f3b4b4dSDmitry Karpeev if (osm->N != PETSC_DETERMINE) { 3558f3b4b4dSDmitry Karpeev /* No local subdomains given, but the desired number of total subdomains is known, so construct them accordingly. */ 3568f3b4b4dSDmitry Karpeev ierr = PCGASMCreateSubdomains(pc->pmat,osm->N,&osm->n,&osm->iis);CHKERRQ(ierr); 3578f3b4b4dSDmitry Karpeev } else if (osm->dm_subdomains && pc->dm) { 3588f3b4b4dSDmitry Karpeev /* try pc->dm next, if allowed */ 3598f3b4b4dSDmitry Karpeev PetscInt d; 360a35b7b57SDmitry Karpeev IS *inner_subdomain_is, *outer_subdomain_is; 361a35b7b57SDmitry Karpeev ierr = DMCreateDomainDecomposition(pc->dm, &num_subdomains, &subdomain_names, &inner_subdomain_is, &outer_subdomain_is, &subdomain_dm);CHKERRQ(ierr); 362a35b7b57SDmitry Karpeev if (num_subdomains) { 363a35b7b57SDmitry Karpeev ierr = PCGASMSetSubdomains(pc, num_subdomains, inner_subdomain_is, outer_subdomain_is);CHKERRQ(ierr); 36469ca1f37SDmitry Karpeev } 365a35b7b57SDmitry Karpeev for (d = 0; d < num_subdomains; ++d) { 366a35b7b57SDmitry Karpeev if (inner_subdomain_is) {ierr = ISDestroy(&inner_subdomain_is[d]);CHKERRQ(ierr);} 367a35b7b57SDmitry Karpeev if (outer_subdomain_is) {ierr = ISDestroy(&outer_subdomain_is[d]);CHKERRQ(ierr);} 368fdc48646SDmitry Karpeev } 369a35b7b57SDmitry Karpeev ierr = PetscFree(inner_subdomain_is);CHKERRQ(ierr); 370a35b7b57SDmitry Karpeev ierr = PetscFree(outer_subdomain_is);CHKERRQ(ierr); 3718f3b4b4dSDmitry Karpeev } else { 3728f3b4b4dSDmitry Karpeev /* still no subdomains; use one per processor */ 37344fe18e5SDmitry Karpeev osm->nmax = osm->n = 1; 374ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr); 375f746d493SDmitry Karpeev osm->N = size; 3768f3b4b4dSDmitry Karpeev ierr = PCGASMCreateLocalSubdomains(pc->pmat,osm->n,&osm->iis);CHKERRQ(ierr); 377fdc48646SDmitry Karpeev } 37806b43e7eSDmitry Karpeev } 379a35b7b57SDmitry Karpeev if (!osm->iis) { 380a35b7b57SDmitry Karpeev /* 3818f3b4b4dSDmitry Karpeev osm->n was set in PCGASMSetSubdomains(), but the actual subdomains have not been supplied. 3828f3b4b4dSDmitry Karpeev We create the requisite number of local inner subdomains and then expand them into 3838f3b4b4dSDmitry Karpeev out subdomains, if necessary. 384a35b7b57SDmitry Karpeev */ 3858f3b4b4dSDmitry Karpeev ierr = PCGASMCreateLocalSubdomains(pc->pmat,osm->n,&osm->iis);CHKERRQ(ierr); 386a35b7b57SDmitry Karpeev } 3878f3b4b4dSDmitry Karpeev if (!osm->ois) { 3888f3b4b4dSDmitry Karpeev /* 3898f3b4b4dSDmitry Karpeev Initially make outer subdomains the same as inner subdomains. If nonzero additional overlap 3908f3b4b4dSDmitry Karpeev has been requested, copy the inner subdomains over so they can be modified. 3918f3b4b4dSDmitry Karpeev */ 3928f3b4b4dSDmitry Karpeev ierr = PetscMalloc1(osm->n,&osm->ois);CHKERRQ(ierr); 3938f3b4b4dSDmitry Karpeev for (i=0; i<osm->n; ++i) { 394ea91fabdSFande Kong if (osm->overlap > 0 && osm->N>1) { /* With positive overlap, osm->iis[i] will be modified */ 3958f3b4b4dSDmitry Karpeev ierr = ISDuplicate(osm->iis[i],(osm->ois)+i);CHKERRQ(ierr); 3968f3b4b4dSDmitry Karpeev ierr = ISCopy(osm->iis[i],osm->ois[i]);CHKERRQ(ierr); 3978f3b4b4dSDmitry Karpeev } else { 3988f3b4b4dSDmitry Karpeev ierr = PetscObjectReference((PetscObject)((osm->iis)[i]));CHKERRQ(ierr); 3998f3b4b4dSDmitry Karpeev osm->ois[i] = osm->iis[i]; 4008f3b4b4dSDmitry Karpeev } 4018f3b4b4dSDmitry Karpeev } 402ea91fabdSFande Kong if (osm->overlap>0 && osm->N>1) { 4038f3b4b4dSDmitry Karpeev /* Extend the "overlapping" regions by a number of steps */ 404d21f2a47SFande Kong ierr = MatIncreaseOverlapSplit(pc->pmat,osm->n,osm->ois,osm->overlap);CHKERRQ(ierr); 4058f3b4b4dSDmitry Karpeev } 406b1a0a8a3SJed Brown } 4076a4f0f73SDmitry Karpeev 4088f3b4b4dSDmitry Karpeev /* Now the subdomains are defined. Determine their global and max local numbers, if necessary. */ 4098f3b4b4dSDmitry Karpeev if (osm->nmax == PETSC_DETERMINE) { 4108f3b4b4dSDmitry Karpeev PetscMPIInt inwork,outwork; 4118f3b4b4dSDmitry Karpeev /* determine global number of subdomains and the max number of local subdomains */ 4128f3b4b4dSDmitry Karpeev inwork = osm->n; 413b2566f29SBarry Smith ierr = MPIU_Allreduce(&inwork,&outwork,1,MPI_INT,MPI_MAX,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr); 4148f3b4b4dSDmitry Karpeev osm->nmax = outwork; 4158f3b4b4dSDmitry Karpeev } 4168f3b4b4dSDmitry Karpeev if (osm->N == PETSC_DETERMINE) { 4178f3b4b4dSDmitry Karpeev /* Determine the number of globally-distinct subdomains and compute a global numbering for them. */ 4188f3b4b4dSDmitry Karpeev ierr = PetscObjectsListGetGlobalNumbering(PetscObjectComm((PetscObject)pc),osm->n,(PetscObject*)osm->ois,&osm->N,NULL);CHKERRQ(ierr); 4198f3b4b4dSDmitry Karpeev } 4208f3b4b4dSDmitry Karpeev 421b1a0a8a3SJed Brown 422b1a0a8a3SJed Brown if (osm->sort_indices) { 423f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 4246a4f0f73SDmitry Karpeev ierr = ISSort(osm->ois[i]);CHKERRQ(ierr); 4256a4f0f73SDmitry Karpeev ierr = ISSort(osm->iis[i]);CHKERRQ(ierr); 426b1a0a8a3SJed Brown } 427b1a0a8a3SJed Brown } 4288f3b4b4dSDmitry Karpeev ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 4298f3b4b4dSDmitry Karpeev ierr = PCGASMPrintSubdomains(pc);CHKERRQ(ierr); 4308f3b4b4dSDmitry Karpeev 4316a4f0f73SDmitry Karpeev /* 4326a4f0f73SDmitry Karpeev Merge the ISs, create merged vectors and restrictions. 4336a4f0f73SDmitry Karpeev */ 4346a4f0f73SDmitry Karpeev /* Merge outer subdomain ISs and construct a restriction onto the disjoint union of local outer subdomains. */ 4356a4f0f73SDmitry Karpeev on = 0; 436f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 4376a4f0f73SDmitry Karpeev ierr = ISGetLocalSize(osm->ois[i],&oni);CHKERRQ(ierr); 4386a4f0f73SDmitry Karpeev on += oni; 439f746d493SDmitry Karpeev } 440785e854fSJed Brown ierr = PetscMalloc1(on, &oidx);CHKERRQ(ierr); 4416a4f0f73SDmitry Karpeev on = 0; 442930d09c1SFande Kong /* Merge local indices together */ 443f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 4446a4f0f73SDmitry Karpeev ierr = ISGetLocalSize(osm->ois[i],&oni);CHKERRQ(ierr); 4456a4f0f73SDmitry Karpeev ierr = ISGetIndices(osm->ois[i],&oidxi);CHKERRQ(ierr); 4466a4f0f73SDmitry Karpeev ierr = PetscMemcpy(oidx+on,oidxi,sizeof(PetscInt)*oni);CHKERRQ(ierr); 4476a4f0f73SDmitry Karpeev ierr = ISRestoreIndices(osm->ois[i],&oidxi);CHKERRQ(ierr); 4486a4f0f73SDmitry Karpeev on += oni; 449f746d493SDmitry Karpeev } 4506a4f0f73SDmitry Karpeev ierr = ISCreateGeneral(((PetscObject)(pc))->comm,on,oidx,PETSC_OWN_POINTER,&gois);CHKERRQ(ierr); 451ea91fabdSFande Kong nTotalInnerIndices = 0; 452ea91fabdSFande Kong for(i=0; i<osm->n; i++){ 453ea91fabdSFande Kong ierr = ISGetLocalSize(osm->iis[i],&nInnerIndices);CHKERRQ(ierr); 454ea91fabdSFande Kong nTotalInnerIndices += nInnerIndices; 455ea91fabdSFande Kong } 456ea91fabdSFande Kong ierr = VecCreateMPI(((PetscObject)(pc))->comm,nTotalInnerIndices,PETSC_DETERMINE,&x);CHKERRQ(ierr); 457ea91fabdSFande Kong ierr = VecDuplicate(x,&y);CHKERRQ(ierr); 458ea91fabdSFande Kong 459ce94432eSBarry Smith ierr = VecCreateMPI(PetscObjectComm((PetscObject)pc),on,PETSC_DECIDE,&osm->gx);CHKERRQ(ierr); 460f746d493SDmitry Karpeev ierr = VecDuplicate(osm->gx,&osm->gy);CHKERRQ(ierr); 461930d09c1SFande Kong ierr = VecGetOwnershipRange(osm->gx, &gostart, NULL);CHKERRQ(ierr); 462930d09c1SFande Kong ierr = ISCreateStride(PetscObjectComm((PetscObject)pc),on,gostart,1, &goid);CHKERRQ(ierr); 463930d09c1SFande Kong /* gois might indices not on local */ 464*9448b7f1SJunchao Zhang ierr = VecScatterCreate(x,gois,osm->gx,goid, &(osm->gorestriction));CHKERRQ(ierr); 465f771a274SFande Kong ierr = PetscMalloc1(osm->n,&numbering);CHKERRQ(ierr); 466f771a274SFande Kong ierr = PetscObjectsListGetGlobalNumbering(PetscObjectComm((PetscObject)pc),osm->n,(PetscObject*)osm->ois,NULL,numbering);CHKERRQ(ierr); 4676a4f0f73SDmitry Karpeev ierr = VecDestroy(&x);CHKERRQ(ierr); 4687105d80bSDmitry Karpeev ierr = ISDestroy(&gois);CHKERRQ(ierr); 4692fa5cd67SKarl Rupp 4706a4f0f73SDmitry Karpeev /* Merge inner subdomain ISs and construct a restriction onto the disjoint union of local inner subdomains. */ 4712fa5cd67SKarl Rupp { 4722fa5cd67SKarl Rupp PetscInt ini; /* Number of indices the i-th a local inner subdomain. */ 4736a4f0f73SDmitry Karpeev PetscInt in; /* Number of indices in the disjoint uniont of local inner subdomains. */ 4746a4f0f73SDmitry Karpeev PetscInt *iidx; /* Global indices in the merged local inner subdomain. */ 4756a4f0f73SDmitry Karpeev PetscInt *ioidx; /* Global indices of the disjoint union of inner subdomains within the disjoint union of outer subdomains. */ 4766a4f0f73SDmitry Karpeev IS giis; /* IS for the disjoint union of inner subdomains. */ 4776a4f0f73SDmitry Karpeev IS giois; /* IS for the disjoint union of inner subdomains within the disjoint union of outer subdomains. */ 478f771a274SFande Kong PetscScalar *array; 479f771a274SFande Kong const PetscInt *indices; 480f771a274SFande Kong PetscInt k; 4816a4f0f73SDmitry Karpeev on = 0; 482f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 4836a4f0f73SDmitry Karpeev ierr = ISGetLocalSize(osm->ois[i],&oni);CHKERRQ(ierr); 4846a4f0f73SDmitry Karpeev on += oni; 485f746d493SDmitry Karpeev } 486f771a274SFande Kong ierr = PetscMalloc1(on, &iidx);CHKERRQ(ierr); 487f771a274SFande Kong ierr = PetscMalloc1(on, &ioidx);CHKERRQ(ierr); 488f771a274SFande Kong ierr = VecGetArray(y,&array);CHKERRQ(ierr); 489e12b4729SFande Kong /* set communicator id to determine where overlap is */ 490f771a274SFande Kong in = 0; 491f771a274SFande Kong for (i=0; i<osm->n; i++) { 492f771a274SFande Kong ierr = ISGetLocalSize(osm->iis[i],&ini);CHKERRQ(ierr); 493f771a274SFande Kong for (k = 0; k < ini; ++k){ 494f771a274SFande Kong array[in+k] = numbering[i]; 495f771a274SFande Kong } 496f771a274SFande Kong in += ini; 497f771a274SFande Kong } 498f771a274SFande Kong ierr = VecRestoreArray(y,&array);CHKERRQ(ierr); 499f771a274SFande Kong ierr = VecScatterBegin(osm->gorestriction,y,osm->gy,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 500f771a274SFande Kong ierr = VecScatterEnd(osm->gorestriction,y,osm->gy,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 501f771a274SFande Kong ierr = VecGetOwnershipRange(osm->gy,&gostart, NULL);CHKERRQ(ierr); 502f771a274SFande Kong ierr = VecGetArray(osm->gy,&array);CHKERRQ(ierr); 503f771a274SFande Kong on = 0; 504f771a274SFande Kong in = 0; 505f771a274SFande Kong for (i=0; i<osm->n; i++) { 506f771a274SFande Kong ierr = ISGetLocalSize(osm->ois[i],&oni);CHKERRQ(ierr); 507f771a274SFande Kong ierr = ISGetIndices(osm->ois[i],&indices);CHKERRQ(ierr); 508f771a274SFande Kong for (k=0; k<oni; k++) { 509e12b4729SFande Kong /* skip overlapping indices to get inner domain */ 51043081b6cSDmitry Karpeev if(PetscRealPart(array[on+k]) != numbering[i]) continue; 511f771a274SFande Kong iidx[in] = indices[k]; 512f771a274SFande Kong ioidx[in++] = gostart+on+k; 513f771a274SFande Kong } 514f771a274SFande Kong ierr = ISRestoreIndices(osm->ois[i], &indices);CHKERRQ(ierr); 515f771a274SFande Kong on += oni; 516f771a274SFande Kong } 517f771a274SFande Kong ierr = VecRestoreArray(osm->gy,&array);CHKERRQ(ierr); 518ce94432eSBarry Smith ierr = ISCreateGeneral(PetscObjectComm((PetscObject)pc),in,iidx,PETSC_OWN_POINTER,&giis);CHKERRQ(ierr); 519ce94432eSBarry Smith ierr = ISCreateGeneral(PetscObjectComm((PetscObject)pc),in,ioidx,PETSC_OWN_POINTER,&giois);CHKERRQ(ierr); 520*9448b7f1SJunchao Zhang ierr = VecScatterCreate(y,giis,osm->gy,giois,&osm->girestriction);CHKERRQ(ierr); 5216a4f0f73SDmitry Karpeev ierr = VecDestroy(&y);CHKERRQ(ierr); 5226a4f0f73SDmitry Karpeev ierr = ISDestroy(&giis);CHKERRQ(ierr); 5236a4f0f73SDmitry Karpeev ierr = ISDestroy(&giois);CHKERRQ(ierr); 524b1a0a8a3SJed Brown } 5256a4f0f73SDmitry Karpeev ierr = ISDestroy(&goid);CHKERRQ(ierr); 526f771a274SFande Kong ierr = PetscFree(numbering);CHKERRQ(ierr); 5272fa5cd67SKarl Rupp 528f746d493SDmitry Karpeev /* Create the subdomain work vectors. */ 529785e854fSJed Brown ierr = PetscMalloc1(osm->n,&osm->x);CHKERRQ(ierr); 530785e854fSJed Brown ierr = PetscMalloc1(osm->n,&osm->y);CHKERRQ(ierr); 531f746d493SDmitry Karpeev ierr = VecGetArray(osm->gx, &gxarray);CHKERRQ(ierr); 532f746d493SDmitry Karpeev ierr = VecGetArray(osm->gy, &gyarray);CHKERRQ(ierr); 5336a4f0f73SDmitry Karpeev for (i=0, on=0; i<osm->n; ++i, on += oni) { 5346a4f0f73SDmitry Karpeev PetscInt oNi; 5356a4f0f73SDmitry Karpeev ierr = ISGetLocalSize(osm->ois[i],&oni);CHKERRQ(ierr); 536930d09c1SFande Kong /* on a sub communicator */ 5376a4f0f73SDmitry Karpeev ierr = ISGetSize(osm->ois[i],&oNi);CHKERRQ(ierr); 5386a4f0f73SDmitry Karpeev ierr = VecCreateMPIWithArray(((PetscObject)(osm->ois[i]))->comm,1,oni,oNi,gxarray+on,&osm->x[i]);CHKERRQ(ierr); 5396a4f0f73SDmitry Karpeev ierr = VecCreateMPIWithArray(((PetscObject)(osm->ois[i]))->comm,1,oni,oNi,gyarray+on,&osm->y[i]);CHKERRQ(ierr); 540b1a0a8a3SJed Brown } 541f746d493SDmitry Karpeev ierr = VecRestoreArray(osm->gx, &gxarray);CHKERRQ(ierr); 542f746d493SDmitry Karpeev ierr = VecRestoreArray(osm->gy, &gyarray);CHKERRQ(ierr); 5438f3b4b4dSDmitry Karpeev /* Create the subdomain solvers */ 544785e854fSJed Brown ierr = PetscMalloc1(osm->n,&osm->ksp);CHKERRQ(ierr); 5458f3b4b4dSDmitry Karpeev for (i=0; i<osm->n; i++) { 5468f3b4b4dSDmitry Karpeev char subprefix[PETSC_MAX_PATH_LEN+1]; 5476a4f0f73SDmitry Karpeev ierr = KSPCreate(((PetscObject)(osm->ois[i]))->comm,&ksp);CHKERRQ(ierr); 54866b14c46SBarry Smith ierr = KSPSetErrorIfNotConverged(ksp,pc->erroriffailure);CHKERRQ(ierr); 5493bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)ksp);CHKERRQ(ierr); 550b1a0a8a3SJed Brown ierr = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,1);CHKERRQ(ierr); 551b1a0a8a3SJed Brown ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr); 5528f3b4b4dSDmitry Karpeev ierr = KSPGetPC(ksp,&subpc);CHKERRQ(ierr); /* Why do we need this here? */ 5538f3b4b4dSDmitry Karpeev if (subdomain_dm) { 5548f3b4b4dSDmitry Karpeev ierr = KSPSetDM(ksp,subdomain_dm[i]);CHKERRQ(ierr); 5558f3b4b4dSDmitry Karpeev ierr = DMDestroy(subdomain_dm+i);CHKERRQ(ierr); 5568f3b4b4dSDmitry Karpeev } 557b1a0a8a3SJed Brown ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 558b1a0a8a3SJed Brown ierr = KSPSetOptionsPrefix(ksp,prefix);CHKERRQ(ierr); 5598f3b4b4dSDmitry Karpeev if (subdomain_names && subdomain_names[i]) { 5608f3b4b4dSDmitry Karpeev ierr = PetscSNPrintf(subprefix,PETSC_MAX_PATH_LEN,"sub_%s_",subdomain_names[i]);CHKERRQ(ierr); 5618f3b4b4dSDmitry Karpeev ierr = KSPAppendOptionsPrefix(ksp,subprefix);CHKERRQ(ierr); 5628f3b4b4dSDmitry Karpeev ierr = PetscFree(subdomain_names[i]);CHKERRQ(ierr); 5638f3b4b4dSDmitry Karpeev } 564b1a0a8a3SJed Brown ierr = KSPAppendOptionsPrefix(ksp,"sub_");CHKERRQ(ierr); 565b1a0a8a3SJed Brown osm->ksp[i] = ksp; 566b1a0a8a3SJed Brown } 5678f3b4b4dSDmitry Karpeev ierr = PetscFree(subdomain_dm);CHKERRQ(ierr); 5688f3b4b4dSDmitry Karpeev ierr = PetscFree(subdomain_names);CHKERRQ(ierr); 569b1a0a8a3SJed Brown scall = MAT_INITIAL_MATRIX; 570b1a0a8a3SJed Brown 5718f3b4b4dSDmitry Karpeev } else { /* if (pc->setupcalled) */ 572b1a0a8a3SJed Brown /* 5738f3b4b4dSDmitry Karpeev Destroy the submatrices from the previous iteration 574b1a0a8a3SJed Brown */ 575b1a0a8a3SJed Brown if (pc->flag == DIFFERENT_NONZERO_PATTERN) { 576f746d493SDmitry Karpeev ierr = MatDestroyMatrices(osm->n,&osm->pmat);CHKERRQ(ierr); 577b1a0a8a3SJed Brown scall = MAT_INITIAL_MATRIX; 578b1a0a8a3SJed Brown } 579ea91fabdSFande Kong if(osm->permutationIS){ 5807dae84e0SHong Zhang ierr = MatCreateSubMatrix(pc->pmat,osm->permutationIS,osm->permutationIS,scall,&osm->permutationP);CHKERRQ(ierr); 581ea91fabdSFande Kong ierr = PetscObjectReference((PetscObject)osm->permutationP);CHKERRQ(ierr); 582ea91fabdSFande Kong osm->pcmat = pc->pmat; 583ea91fabdSFande Kong pc->pmat = osm->permutationP; 584b1a0a8a3SJed Brown } 585b1a0a8a3SJed Brown 586ea91fabdSFande Kong } 587ea91fabdSFande Kong 588ea91fabdSFande Kong 589b1a0a8a3SJed Brown /* 590f746d493SDmitry Karpeev Extract out the submatrices. 591b1a0a8a3SJed Brown */ 59281a5c4aaSDmitry Karpeev if (size > 1) { 5937dae84e0SHong Zhang ierr = MatCreateSubMatricesMPI(pc->pmat,osm->n,osm->ois,osm->ois,scall,&osm->pmat);CHKERRQ(ierr); 5942fa5cd67SKarl Rupp } else { 5957dae84e0SHong Zhang ierr = MatCreateSubMatrices(pc->pmat,osm->n,osm->ois,osm->ois,scall,&osm->pmat);CHKERRQ(ierr); 59681a5c4aaSDmitry Karpeev } 597b1a0a8a3SJed Brown if (scall == MAT_INITIAL_MATRIX) { 598b1a0a8a3SJed Brown ierr = PetscObjectGetOptionsPrefix((PetscObject)pc->pmat,&pprefix);CHKERRQ(ierr); 599f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 6003bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)osm->pmat[i]);CHKERRQ(ierr); 601b1a0a8a3SJed Brown ierr = PetscObjectSetOptionsPrefix((PetscObject)osm->pmat[i],pprefix);CHKERRQ(ierr); 602b1a0a8a3SJed Brown } 603b1a0a8a3SJed Brown } 604b1a0a8a3SJed Brown 605b1a0a8a3SJed Brown /* Return control to the user so that the submatrices can be modified (e.g., to apply 606b1a0a8a3SJed Brown different boundary conditions for the submatrices than for the global problem) */ 6076a4f0f73SDmitry Karpeev ierr = PCModifySubMatrices(pc,osm->n,osm->ois,osm->ois,osm->pmat,pc->modifysubmatricesP);CHKERRQ(ierr); 608b1a0a8a3SJed Brown 609b1a0a8a3SJed Brown /* 6106a4f0f73SDmitry Karpeev Loop over submatrices putting them into local ksps 611b1a0a8a3SJed Brown */ 612f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 61323ee1639SBarry Smith ierr = KSPSetOperators(osm->ksp[i],osm->pmat[i],osm->pmat[i]);CHKERRQ(ierr); 614b1a0a8a3SJed Brown if (!pc->setupcalled) { 615b1a0a8a3SJed Brown ierr = KSPSetFromOptions(osm->ksp[i]);CHKERRQ(ierr); 616b1a0a8a3SJed Brown } 617b1a0a8a3SJed Brown } 618ea91fabdSFande Kong if(osm->pcmat){ 619ea91fabdSFande Kong ierr = MatDestroy(&pc->pmat);CHKERRQ(ierr); 620ea91fabdSFande Kong pc->pmat = osm->pcmat; 621ea91fabdSFande Kong osm->pcmat = 0; 622ea91fabdSFande Kong } 623b1a0a8a3SJed Brown PetscFunctionReturn(0); 624b1a0a8a3SJed Brown } 625b1a0a8a3SJed Brown 626f746d493SDmitry Karpeev static PetscErrorCode PCSetUpOnBlocks_GASM(PC pc) 627b1a0a8a3SJed Brown { 628f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 629b1a0a8a3SJed Brown PetscErrorCode ierr; 630b1a0a8a3SJed Brown PetscInt i; 631b1a0a8a3SJed Brown 632b1a0a8a3SJed Brown PetscFunctionBegin; 633f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 634b1a0a8a3SJed Brown ierr = KSPSetUp(osm->ksp[i]);CHKERRQ(ierr); 635b1a0a8a3SJed Brown } 636b1a0a8a3SJed Brown PetscFunctionReturn(0); 637b1a0a8a3SJed Brown } 638b1a0a8a3SJed Brown 639ea91fabdSFande Kong static PetscErrorCode PCApply_GASM(PC pc,Vec xin,Vec yout) 640b1a0a8a3SJed Brown { 641f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 642b1a0a8a3SJed Brown PetscErrorCode ierr; 643f746d493SDmitry Karpeev PetscInt i; 644ea91fabdSFande Kong Vec x,y; 645b1a0a8a3SJed Brown ScatterMode forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE; 646b1a0a8a3SJed Brown 647b1a0a8a3SJed Brown PetscFunctionBegin; 648ea91fabdSFande Kong if(osm->pctoouter){ 649ea91fabdSFande Kong ierr = VecScatterBegin(osm->pctoouter,xin,osm->pcx,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 650ea91fabdSFande Kong ierr = VecScatterEnd(osm->pctoouter,xin,osm->pcx,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 651ea91fabdSFande Kong x = osm->pcx; 652ea91fabdSFande Kong y = osm->pcy; 653ea91fabdSFande Kong }else{ 654ea91fabdSFande Kong x = xin; 655ea91fabdSFande Kong y = yout; 656ea91fabdSFande Kong } 657b1a0a8a3SJed Brown /* 6586a4f0f73SDmitry Karpeev Support for limiting the restriction or interpolation only to the inner 659b1a0a8a3SJed Brown subdomain values (leaving the other values 0). 660b1a0a8a3SJed Brown */ 661f746d493SDmitry Karpeev if (!(osm->type & PC_GASM_RESTRICT)) { 662b1a0a8a3SJed Brown /* have to zero the work RHS since scatter may leave some slots empty */ 663f746d493SDmitry Karpeev ierr = VecZeroEntries(osm->gx);CHKERRQ(ierr); 6646a4f0f73SDmitry Karpeev ierr = VecScatterBegin(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 6652fa5cd67SKarl Rupp } else { 6666a4f0f73SDmitry Karpeev ierr = VecScatterBegin(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 667b1a0a8a3SJed Brown } 6686a4f0f73SDmitry Karpeev ierr = VecZeroEntries(osm->gy);CHKERRQ(ierr); 6696a4f0f73SDmitry Karpeev if (!(osm->type & PC_GASM_RESTRICT)) { 6706a4f0f73SDmitry Karpeev ierr = VecScatterEnd(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 6712fa5cd67SKarl Rupp } else { 6726a4f0f73SDmitry Karpeev ierr = VecScatterEnd(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 6736a4f0f73SDmitry Karpeev } 67486b96d36SDmitry Karpeev /* do the subdomain solves */ 67586b96d36SDmitry Karpeev for (i=0; i<osm->n; ++i) { 676b1a0a8a3SJed Brown ierr = KSPSolve(osm->ksp[i],osm->x[i],osm->y[i]);CHKERRQ(ierr); 677c0decd05SBarry Smith ierr = KSPCheckSolve(osm->ksp[i],pc,osm->y[i]);CHKERRQ(ierr); 678b1a0a8a3SJed Brown } 679930d09c1SFande Kong /* Do we need to zero y ?? */ 6806a4f0f73SDmitry Karpeev ierr = VecZeroEntries(y);CHKERRQ(ierr); 6816a4f0f73SDmitry Karpeev if (!(osm->type & PC_GASM_INTERPOLATE)) { 6826a4f0f73SDmitry Karpeev ierr = VecScatterBegin(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); 683ea91fabdSFande Kong ierr = VecScatterEnd(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); 6842fa5cd67SKarl Rupp } else { 6856a4f0f73SDmitry Karpeev ierr = VecScatterBegin(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); 686ea91fabdSFande Kong ierr = VecScatterEnd(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); 6876a4f0f73SDmitry Karpeev } 688ea91fabdSFande Kong if(osm->pctoouter){ 689ea91fabdSFande Kong ierr = VecScatterBegin(osm->pctoouter,y,yout,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 690ea91fabdSFande Kong ierr = VecScatterEnd(osm->pctoouter,y,yout,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 691ea91fabdSFande Kong } 692ea91fabdSFande Kong PetscFunctionReturn(0); 693b1a0a8a3SJed Brown } 694b1a0a8a3SJed Brown 695ea91fabdSFande Kong static PetscErrorCode PCApplyTranspose_GASM(PC pc,Vec xin,Vec yout) 696b1a0a8a3SJed Brown { 697f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 698b1a0a8a3SJed Brown PetscErrorCode ierr; 699f746d493SDmitry Karpeev PetscInt i; 700ea91fabdSFande Kong Vec x,y; 701b1a0a8a3SJed Brown ScatterMode forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE; 702b1a0a8a3SJed Brown 703b1a0a8a3SJed Brown PetscFunctionBegin; 704ea91fabdSFande Kong if(osm->pctoouter){ 705ea91fabdSFande Kong ierr = VecScatterBegin(osm->pctoouter,xin,osm->pcx,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 706ea91fabdSFande Kong ierr = VecScatterEnd(osm->pctoouter,xin,osm->pcx,INSERT_VALUES,SCATTER_REVERSE);CHKERRQ(ierr); 707ea91fabdSFande Kong x = osm->pcx; 708ea91fabdSFande Kong y = osm->pcy; 709ea91fabdSFande Kong }else{ 710ea91fabdSFande Kong x = xin; 711ea91fabdSFande Kong y = yout; 712ea91fabdSFande Kong } 713b1a0a8a3SJed Brown /* 714b1a0a8a3SJed Brown Support for limiting the restriction or interpolation to only local 715b1a0a8a3SJed Brown subdomain values (leaving the other values 0). 716b1a0a8a3SJed Brown 717f746d493SDmitry Karpeev Note: these are reversed from the PCApply_GASM() because we are applying the 718b1a0a8a3SJed Brown transpose of the three terms 719b1a0a8a3SJed Brown */ 720f746d493SDmitry Karpeev if (!(osm->type & PC_GASM_INTERPOLATE)) { 721b1a0a8a3SJed Brown /* have to zero the work RHS since scatter may leave some slots empty */ 722f746d493SDmitry Karpeev ierr = VecZeroEntries(osm->gx);CHKERRQ(ierr); 7236a4f0f73SDmitry Karpeev ierr = VecScatterBegin(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 7242fa5cd67SKarl Rupp } else { 7256a4f0f73SDmitry Karpeev ierr = VecScatterBegin(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 726b1a0a8a3SJed Brown } 7276a4f0f73SDmitry Karpeev ierr = VecZeroEntries(osm->gy);CHKERRQ(ierr); 7286a4f0f73SDmitry Karpeev if (!(osm->type & PC_GASM_INTERPOLATE)) { 7296a4f0f73SDmitry Karpeev ierr = VecScatterEnd(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 7302fa5cd67SKarl Rupp } else { 7316a4f0f73SDmitry Karpeev ierr = VecScatterEnd(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 7326a4f0f73SDmitry Karpeev } 733b1a0a8a3SJed Brown /* do the local solves */ 734ab3e923fSDmitry 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. */ 735b1a0a8a3SJed Brown ierr = KSPSolveTranspose(osm->ksp[i],osm->x[i],osm->y[i]);CHKERRQ(ierr); 736c0decd05SBarry Smith ierr = KSPCheckSolve(osm->ksp[i],pc,osm->y[i]);CHKERRQ(ierr); 737b1a0a8a3SJed Brown } 7386a4f0f73SDmitry Karpeev ierr = VecZeroEntries(y);CHKERRQ(ierr); 7396a4f0f73SDmitry Karpeev if (!(osm->type & PC_GASM_RESTRICT)) { 7406a4f0f73SDmitry Karpeev ierr = VecScatterBegin(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); 7416a4f0f73SDmitry Karpeev ierr = VecScatterEnd(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); 7422fa5cd67SKarl Rupp } else { 7436a4f0f73SDmitry Karpeev ierr = VecScatterBegin(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); 7446a4f0f73SDmitry Karpeev ierr = VecScatterEnd(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); 7456a4f0f73SDmitry Karpeev } 746ea91fabdSFande Kong if(osm->pctoouter){ 747ea91fabdSFande Kong ierr = VecScatterBegin(osm->pctoouter,y,yout,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 748ea91fabdSFande Kong ierr = VecScatterEnd(osm->pctoouter,y,yout,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 749ea91fabdSFande Kong } 750b1a0a8a3SJed Brown PetscFunctionReturn(0); 751b1a0a8a3SJed Brown } 752b1a0a8a3SJed Brown 753a06653b4SBarry Smith static PetscErrorCode PCReset_GASM(PC pc) 754b1a0a8a3SJed Brown { 755f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 756b1a0a8a3SJed Brown PetscErrorCode ierr; 757b1a0a8a3SJed Brown PetscInt i; 758b1a0a8a3SJed Brown 759b1a0a8a3SJed Brown PetscFunctionBegin; 760b1a0a8a3SJed Brown if (osm->ksp) { 761f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 762a06653b4SBarry Smith ierr = KSPReset(osm->ksp[i]);CHKERRQ(ierr); 763b1a0a8a3SJed Brown } 764b1a0a8a3SJed Brown } 765b1a0a8a3SJed Brown if (osm->pmat) { 766f746d493SDmitry Karpeev if (osm->n > 0) { 767df750dc8SHong Zhang PetscMPIInt size; 768df750dc8SHong Zhang ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr); 769df750dc8SHong Zhang if (size > 1) { 7707dae84e0SHong Zhang /* osm->pmat is created by MatCreateSubMatricesMPI(), cannot use MatDestroySubMatrices() */ 771ab3e923fSDmitry Karpeev ierr = MatDestroyMatrices(osm->n,&osm->pmat);CHKERRQ(ierr); 772df750dc8SHong Zhang } else { 773df750dc8SHong Zhang ierr = MatDestroySubMatrices(osm->n,&osm->pmat);CHKERRQ(ierr); 774df750dc8SHong Zhang } 775b1a0a8a3SJed Brown } 776b1a0a8a3SJed Brown } 777d34fcf5fSBarry Smith if (osm->x) { 778f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 779fcfd50ebSBarry Smith ierr = VecDestroy(&osm->x[i]);CHKERRQ(ierr); 780fcfd50ebSBarry Smith ierr = VecDestroy(&osm->y[i]);CHKERRQ(ierr); 781b1a0a8a3SJed Brown } 782d34fcf5fSBarry Smith } 7836bf464f9SBarry Smith ierr = VecDestroy(&osm->gx);CHKERRQ(ierr); 7846bf464f9SBarry Smith ierr = VecDestroy(&osm->gy);CHKERRQ(ierr); 785ab3e923fSDmitry Karpeev 7866a4f0f73SDmitry Karpeev ierr = VecScatterDestroy(&osm->gorestriction);CHKERRQ(ierr); 7876a4f0f73SDmitry Karpeev ierr = VecScatterDestroy(&osm->girestriction);CHKERRQ(ierr); 7888f3b4b4dSDmitry Karpeev if (!osm->user_subdomains) { 7892c112581SDmitry Karpeev ierr = PCGASMDestroySubdomains(osm->n,&osm->ois,&osm->iis);CHKERRQ(ierr); 7908f3b4b4dSDmitry Karpeev osm->N = PETSC_DETERMINE; 7918f3b4b4dSDmitry Karpeev osm->nmax = PETSC_DETERMINE; 7928f3b4b4dSDmitry Karpeev } 793ea91fabdSFande Kong if(osm->pctoouter){ 794ea91fabdSFande Kong ierr = VecScatterDestroy(&(osm->pctoouter));CHKERRQ(ierr); 795ea91fabdSFande Kong } 796ea91fabdSFande Kong if(osm->permutationIS){ 797ea91fabdSFande Kong ierr = ISDestroy(&(osm->permutationIS));CHKERRQ(ierr); 798ea91fabdSFande Kong } 799ea91fabdSFande Kong if(osm->pcx){ 800ea91fabdSFande Kong ierr = VecDestroy(&(osm->pcx));CHKERRQ(ierr); 801ea91fabdSFande Kong } 802ea91fabdSFande Kong if(osm->pcy){ 803ea91fabdSFande Kong ierr = VecDestroy(&(osm->pcy));CHKERRQ(ierr); 804ea91fabdSFande Kong } 805ea91fabdSFande Kong if(osm->permutationP){ 806ea91fabdSFande Kong ierr = MatDestroy(&(osm->permutationP));CHKERRQ(ierr); 807ea91fabdSFande Kong } 808ea91fabdSFande Kong if(osm->pcmat){ 809ea91fabdSFande Kong ierr = MatDestroy(&osm->pcmat);CHKERRQ(ierr); 810ea91fabdSFande Kong } 811a06653b4SBarry Smith PetscFunctionReturn(0); 812a06653b4SBarry Smith } 813a06653b4SBarry Smith 814a06653b4SBarry Smith static PetscErrorCode PCDestroy_GASM(PC pc) 815a06653b4SBarry Smith { 816a06653b4SBarry Smith PC_GASM *osm = (PC_GASM*)pc->data; 817a06653b4SBarry Smith PetscErrorCode ierr; 818a06653b4SBarry Smith PetscInt i; 819a06653b4SBarry Smith 820a06653b4SBarry Smith PetscFunctionBegin; 821a06653b4SBarry Smith ierr = PCReset_GASM(pc);CHKERRQ(ierr); 822135757f6SDmitry Karpeev /* PCReset will not destroy subdomains, if user_subdomains is true. */ 823135757f6SDmitry Karpeev ierr = PCGASMDestroySubdomains(osm->n,&osm->ois,&osm->iis);CHKERRQ(ierr); 824a06653b4SBarry Smith if (osm->ksp) { 825a06653b4SBarry Smith for (i=0; i<osm->n; i++) { 8266bf464f9SBarry Smith ierr = KSPDestroy(&osm->ksp[i]);CHKERRQ(ierr); 827a06653b4SBarry Smith } 828a06653b4SBarry Smith ierr = PetscFree(osm->ksp);CHKERRQ(ierr); 829a06653b4SBarry Smith } 830a06653b4SBarry Smith ierr = PetscFree(osm->x);CHKERRQ(ierr); 831a06653b4SBarry Smith ierr = PetscFree(osm->y);CHKERRQ(ierr); 832c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 833b1a0a8a3SJed Brown PetscFunctionReturn(0); 834b1a0a8a3SJed Brown } 835b1a0a8a3SJed Brown 8364416b707SBarry Smith static PetscErrorCode PCSetFromOptions_GASM(PetscOptionItems *PetscOptionsObject,PC pc) 837a6dfd86eSKarl Rupp { 838f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 839b1a0a8a3SJed Brown PetscErrorCode ierr; 840b1a0a8a3SJed Brown PetscInt blocks,ovl; 841b1a0a8a3SJed Brown PetscBool symset,flg; 842f746d493SDmitry Karpeev PCGASMType gasmtype; 843b1a0a8a3SJed Brown 844b1a0a8a3SJed Brown PetscFunctionBegin; 845b1a0a8a3SJed Brown /* set the type to symmetric if matrix is symmetric */ 846b1a0a8a3SJed Brown if (!osm->type_set && pc->pmat) { 847b1a0a8a3SJed Brown ierr = MatIsSymmetricKnown(pc->pmat,&symset,&flg);CHKERRQ(ierr); 8482fa5cd67SKarl Rupp if (symset && flg) osm->type = PC_GASM_BASIC; 849b1a0a8a3SJed Brown } 850e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"Generalized additive Schwarz options");CHKERRQ(ierr); 8518f3b4b4dSDmitry 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); 8528f3b4b4dSDmitry Karpeev ierr = PetscOptionsInt("-pc_gasm_total_subdomains","Total number of subdomains across communicator","PCGASMSetTotalSubdomains",osm->N,&blocks,&flg);CHKERRQ(ierr); 85365db9045SDmitry Karpeev if (flg) { 8548f3b4b4dSDmitry Karpeev ierr = PCGASMSetTotalSubdomains(pc,blocks);CHKERRQ(ierr); 85565db9045SDmitry Karpeev } 85606b43e7eSDmitry Karpeev ierr = PetscOptionsInt("-pc_gasm_overlap","Number of overlapping degrees of freedom","PCGASMSetOverlap",osm->overlap,&ovl,&flg);CHKERRQ(ierr); 85765db9045SDmitry Karpeev if (flg) { 85865db9045SDmitry Karpeev ierr = PCGASMSetOverlap(pc,ovl);CHKERRQ(ierr); 859d709ea83SDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 86065db9045SDmitry Karpeev } 861b1a0a8a3SJed Brown flg = PETSC_FALSE; 862f746d493SDmitry Karpeev ierr = PetscOptionsEnum("-pc_gasm_type","Type of restriction/extension","PCGASMSetType",PCGASMTypes,(PetscEnum)osm->type,(PetscEnum*)&gasmtype,&flg);CHKERRQ(ierr); 863f746d493SDmitry Karpeev if (flg) {ierr = PCGASMSetType(pc,gasmtype);CHKERRQ(ierr);} 864ea91fabdSFande Kong ierr = PetscOptionsBool("-pc_gasm_use_hierachical_partitioning","use hierarchical partitioning",NULL,osm->hierarchicalpartitioning,&osm->hierarchicalpartitioning,&flg);CHKERRQ(ierr); 865b1a0a8a3SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 866b1a0a8a3SJed Brown PetscFunctionReturn(0); 867b1a0a8a3SJed Brown } 868b1a0a8a3SJed Brown 869b1a0a8a3SJed Brown /*------------------------------------------------------------------------------------*/ 870b1a0a8a3SJed Brown 8718f3b4b4dSDmitry Karpeev /*@ 8728f3b4b4dSDmitry Karpeev PCGASMSetTotalSubdomains - sets the total number of subdomains to use across the 8738f3b4b4dSDmitry Karpeev communicator. 8748f3b4b4dSDmitry Karpeev Logically collective on pc 8758f3b4b4dSDmitry Karpeev 8768f3b4b4dSDmitry Karpeev Input Parameters: 8778f3b4b4dSDmitry Karpeev + pc - the preconditioner 8788f3b4b4dSDmitry Karpeev - N - total number of subdomains 8798f3b4b4dSDmitry Karpeev 8808f3b4b4dSDmitry Karpeev 8818f3b4b4dSDmitry Karpeev Level: beginner 8828f3b4b4dSDmitry Karpeev 8838f3b4b4dSDmitry Karpeev .keywords: PC, ASM, DM, set, subdomains, additive Schwarz 8848f3b4b4dSDmitry Karpeev 8858f3b4b4dSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMSetOverlap() 8868f3b4b4dSDmitry Karpeev PCGASMCreateSubdomains2D() 8878f3b4b4dSDmitry Karpeev @*/ 8888f3b4b4dSDmitry Karpeev PetscErrorCode PCGASMSetTotalSubdomains(PC pc,PetscInt N) 8898f3b4b4dSDmitry Karpeev { 8908f3b4b4dSDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 8918f3b4b4dSDmitry Karpeev PetscMPIInt size,rank; 8928f3b4b4dSDmitry Karpeev PetscErrorCode ierr; 8938f3b4b4dSDmitry Karpeev 8948f3b4b4dSDmitry Karpeev PetscFunctionBegin; 8958f3b4b4dSDmitry 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); 8968f3b4b4dSDmitry Karpeev if (pc->setupcalled) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetTotalSubdomains() should be called before calling PCSetUp()."); 8978f3b4b4dSDmitry Karpeev 8982c112581SDmitry Karpeev ierr = PCGASMDestroySubdomains(osm->n,&osm->iis,&osm->ois);CHKERRQ(ierr); 8998f3b4b4dSDmitry Karpeev osm->ois = osm->iis = NULL; 9008f3b4b4dSDmitry Karpeev 9018f3b4b4dSDmitry Karpeev ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr); 9028f3b4b4dSDmitry Karpeev ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank);CHKERRQ(ierr); 9038f3b4b4dSDmitry Karpeev osm->N = N; 9048f3b4b4dSDmitry Karpeev osm->n = PETSC_DETERMINE; 9058f3b4b4dSDmitry Karpeev osm->nmax = PETSC_DETERMINE; 9068f3b4b4dSDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 9078f3b4b4dSDmitry Karpeev PetscFunctionReturn(0); 9088f3b4b4dSDmitry Karpeev } 9098f3b4b4dSDmitry Karpeev 910670417b2SFande Kong 911b541e6a4SDmitry Karpeev static PetscErrorCode PCGASMSetSubdomains_GASM(PC pc,PetscInt n,IS iis[],IS ois[]) 912b1a0a8a3SJed Brown { 913f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 914b1a0a8a3SJed Brown PetscErrorCode ierr; 915b1a0a8a3SJed Brown PetscInt i; 916b1a0a8a3SJed Brown 917b1a0a8a3SJed Brown PetscFunctionBegin; 9188f3b4b4dSDmitry Karpeev if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Each process must have 1 or more subdomains, got n = %D",n); 9198f3b4b4dSDmitry Karpeev if (pc->setupcalled) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetSubdomains() should be called before calling PCSetUp()."); 920b1a0a8a3SJed Brown 9212c112581SDmitry Karpeev ierr = PCGASMDestroySubdomains(osm->n,&osm->iis,&osm->ois);CHKERRQ(ierr); 9228f3b4b4dSDmitry Karpeev osm->iis = osm->ois = NULL; 9238f3b4b4dSDmitry Karpeev osm->n = n; 9248f3b4b4dSDmitry Karpeev osm->N = PETSC_DETERMINE; 9258f3b4b4dSDmitry Karpeev osm->nmax = PETSC_DETERMINE; 926a35b7b57SDmitry Karpeev if (ois) { 927785e854fSJed Brown ierr = PetscMalloc1(n,&osm->ois);CHKERRQ(ierr); 9288f3b4b4dSDmitry Karpeev for (i=0; i<n; i++) { 9298f3b4b4dSDmitry Karpeev ierr = PetscObjectReference((PetscObject)ois[i]);CHKERRQ(ierr); 9308f3b4b4dSDmitry Karpeev osm->ois[i] = ois[i]; 9318f3b4b4dSDmitry Karpeev } 9328f3b4b4dSDmitry Karpeev /* 9338f3b4b4dSDmitry Karpeev Since the user set the outer subdomains, even if nontrivial overlap was requested via PCGASMSetOverlap(), 9348f3b4b4dSDmitry Karpeev it will be ignored. To avoid confusion later on (e.g., when viewing the PC), the overlap size is set to -1. 9358f3b4b4dSDmitry Karpeev */ 936b1a0a8a3SJed Brown osm->overlap = -1; 937670417b2SFande Kong /* inner subdomains must be provided */ 938670417b2SFande Kong if (!iis) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"inner indices have to be provided \n"); 939670417b2SFande Kong }/* end if */ 940a35b7b57SDmitry Karpeev if (iis) { 941785e854fSJed Brown ierr = PetscMalloc1(n,&osm->iis);CHKERRQ(ierr); 9428f3b4b4dSDmitry Karpeev for (i=0; i<n; i++) { 9438f3b4b4dSDmitry Karpeev ierr = PetscObjectReference((PetscObject)iis[i]);CHKERRQ(ierr); 9448f3b4b4dSDmitry Karpeev osm->iis[i] = iis[i]; 9458f3b4b4dSDmitry Karpeev } 946a35b7b57SDmitry Karpeev if (!ois) { 947390e1bf2SBarry Smith osm->ois = NULL; 948670417b2SFande Kong /* if user does not provide outer indices, we will create the corresponding outer indices using osm->overlap =1 in PCSetUp_GASM */ 949670417b2SFande Kong } 950670417b2SFande Kong } 951670417b2SFande Kong #if defined(PETSC_USE_DEBUG) 952670417b2SFande Kong { 953670417b2SFande Kong PetscInt j,rstart,rend,*covered,lsize; 954670417b2SFande Kong const PetscInt *indices; 955670417b2SFande Kong /* check if the inner indices cover and only cover the local portion of the preconditioning matrix */ 956670417b2SFande Kong ierr = MatGetOwnershipRange(pc->pmat,&rstart,&rend);CHKERRQ(ierr); 957670417b2SFande Kong ierr = PetscCalloc1(rend-rstart,&covered);CHKERRQ(ierr); 958670417b2SFande Kong /* check if the current processor owns indices from others */ 959a35b7b57SDmitry Karpeev for (i=0; i<n; i++) { 960670417b2SFande Kong ierr = ISGetIndices(osm->iis[i],&indices);CHKERRQ(ierr); 961670417b2SFande Kong ierr = ISGetLocalSize(osm->iis[i],&lsize);CHKERRQ(ierr); 962670417b2SFande Kong for (j=0; j<lsize; j++) { 963670417b2SFande Kong if (indices[j]<rstart || indices[j]>=rend) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"inner subdomains can not own an index %d from other processors", indices[j]); 964670417b2SFande Kong else if (covered[indices[j]-rstart]==1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"inner subdomains can not have an overlapping index %d ",indices[j]); 965670417b2SFande Kong else covered[indices[j]-rstart] = 1; 966a35b7b57SDmitry Karpeev } 967670417b2SFande Kong ierr = ISRestoreIndices(osm->iis[i],&indices);CHKERRQ(ierr); 9688f3b4b4dSDmitry Karpeev } 969670417b2SFande Kong /* check if we miss any indices */ 970670417b2SFande Kong for (i=rstart; i<rend; i++) { 971670417b2SFande Kong if (!covered[i-rstart]) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"local entity %d was not covered by inner subdomains",i); 972a35b7b57SDmitry Karpeev } 973670417b2SFande Kong ierr = PetscFree(covered);CHKERRQ(ierr); 974a35b7b57SDmitry Karpeev } 975670417b2SFande Kong #endif 9768f3b4b4dSDmitry Karpeev if (iis) osm->user_subdomains = PETSC_TRUE; 977b1a0a8a3SJed Brown PetscFunctionReturn(0); 978b1a0a8a3SJed Brown } 979b1a0a8a3SJed Brown 980b1a0a8a3SJed Brown 981f7a08781SBarry Smith static PetscErrorCode PCGASMSetOverlap_GASM(PC pc,PetscInt ovl) 982b1a0a8a3SJed Brown { 983f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 984b1a0a8a3SJed Brown 985b1a0a8a3SJed Brown PetscFunctionBegin; 986ce94432eSBarry Smith if (ovl < 0) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap value requested"); 987ce94432eSBarry Smith if (pc->setupcalled && ovl != osm->overlap) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetOverlap() should be called before PCSetUp()."); 9882fa5cd67SKarl Rupp if (!pc->setupcalled) osm->overlap = ovl; 989b1a0a8a3SJed Brown PetscFunctionReturn(0); 990b1a0a8a3SJed Brown } 991b1a0a8a3SJed Brown 992f7a08781SBarry Smith static PetscErrorCode PCGASMSetType_GASM(PC pc,PCGASMType type) 993b1a0a8a3SJed Brown { 994f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 995b1a0a8a3SJed Brown 996b1a0a8a3SJed Brown PetscFunctionBegin; 997b1a0a8a3SJed Brown osm->type = type; 998b1a0a8a3SJed Brown osm->type_set = PETSC_TRUE; 999b1a0a8a3SJed Brown PetscFunctionReturn(0); 1000b1a0a8a3SJed Brown } 1001b1a0a8a3SJed Brown 1002f7a08781SBarry Smith static PetscErrorCode PCGASMSetSortIndices_GASM(PC pc,PetscBool doSort) 1003b1a0a8a3SJed Brown { 1004f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 1005b1a0a8a3SJed Brown 1006b1a0a8a3SJed Brown PetscFunctionBegin; 1007b1a0a8a3SJed Brown osm->sort_indices = doSort; 1008b1a0a8a3SJed Brown PetscFunctionReturn(0); 1009b1a0a8a3SJed Brown } 1010b1a0a8a3SJed Brown 101144fe18e5SDmitry Karpeev /* 10128f3b4b4dSDmitry Karpeev FIXME: This routine might need to be modified now that multiple ranks per subdomain are allowed. 101344fe18e5SDmitry Karpeev In particular, it would upset the global subdomain number calculation. 101444fe18e5SDmitry Karpeev */ 1015f7a08781SBarry Smith static PetscErrorCode PCGASMGetSubKSP_GASM(PC pc,PetscInt *n,PetscInt *first,KSP **ksp) 1016b1a0a8a3SJed Brown { 1017f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 1018b1a0a8a3SJed Brown PetscErrorCode ierr; 1019b1a0a8a3SJed Brown 1020b1a0a8a3SJed Brown PetscFunctionBegin; 1021ce94432eSBarry 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"); 1022b1a0a8a3SJed Brown 10232fa5cd67SKarl Rupp if (n) *n = osm->n; 1024ab3e923fSDmitry Karpeev if (first) { 1025ce94432eSBarry Smith ierr = MPI_Scan(&osm->n,first,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)pc));CHKERRQ(ierr); 1026ab3e923fSDmitry Karpeev *first -= osm->n; 1027b1a0a8a3SJed Brown } 1028b1a0a8a3SJed Brown if (ksp) { 1029b1a0a8a3SJed Brown /* Assume that local solves are now different; not necessarily 103006b43e7eSDmitry Karpeev true, though! This flag is used only for PCView_GASM() */ 1031b1a0a8a3SJed Brown *ksp = osm->ksp; 10326a4f0f73SDmitry Karpeev osm->same_subdomain_solvers = PETSC_FALSE; 1033b1a0a8a3SJed Brown } 1034b1a0a8a3SJed Brown PetscFunctionReturn(0); 1035ab3e923fSDmitry Karpeev } /* PCGASMGetSubKSP_GASM() */ 1036b1a0a8a3SJed Brown 1037b1a0a8a3SJed Brown /*@C 103806b43e7eSDmitry Karpeev PCGASMSetSubdomains - Sets the subdomains for this processor 103906b43e7eSDmitry Karpeev for the additive Schwarz preconditioner. 1040b1a0a8a3SJed Brown 10418f3b4b4dSDmitry Karpeev Collective on pc 1042b1a0a8a3SJed Brown 1043b1a0a8a3SJed Brown Input Parameters: 10448f3b4b4dSDmitry Karpeev + pc - the preconditioner object 104506b43e7eSDmitry Karpeev . n - the number of subdomains for this processor 10468f3b4b4dSDmitry Karpeev . iis - the index sets that define the inner subdomains (or NULL for PETSc to determine subdomains) 10478f3b4b4dSDmitry 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) 1048b1a0a8a3SJed Brown 1049b1a0a8a3SJed Brown Notes: 105006b43e7eSDmitry Karpeev The IS indices use the parallel, global numbering of the vector entries. 10516a4f0f73SDmitry Karpeev Inner subdomains are those where the correction is applied. 10526a4f0f73SDmitry Karpeev Outer subdomains are those where the residual necessary to obtain the 10536a4f0f73SDmitry Karpeev corrections is obtained (see PCGASMType for the use of inner/outer subdomains). 10546a4f0f73SDmitry Karpeev Both inner and outer subdomains can extend over several processors. 10556a4f0f73SDmitry Karpeev This processor's portion of a subdomain is known as a local subdomain. 10566a4f0f73SDmitry Karpeev 1057670417b2SFande Kong Inner subdomains can not overlap with each other, do not have any entities from remote processors, 1058670417b2SFande Kong and have to cover the entire local subdomain owned by the current processor. The index sets on each 1059670417b2SFande Kong process should be ordered such that the ith local subdomain is connected to the ith remote subdomain 1060670417b2SFande Kong on another MPI process. 1061670417b2SFande Kong 10626a4f0f73SDmitry Karpeev By default the GASM preconditioner uses 1 (local) subdomain per processor. 10636a4f0f73SDmitry Karpeev 1064b1a0a8a3SJed Brown 1065b1a0a8a3SJed Brown Level: advanced 1066b1a0a8a3SJed Brown 106706b43e7eSDmitry Karpeev .keywords: PC, GASM, set, subdomains, additive Schwarz 1068b1a0a8a3SJed Brown 10698f3b4b4dSDmitry Karpeev .seealso: PCGASMSetNumSubdomains(), PCGASMSetOverlap(), PCGASMGetSubKSP(), 107006b43e7eSDmitry Karpeev PCGASMCreateSubdomains2D(), PCGASMGetSubdomains() 1071b1a0a8a3SJed Brown @*/ 10726a4f0f73SDmitry Karpeev PetscErrorCode PCGASMSetSubdomains(PC pc,PetscInt n,IS iis[],IS ois[]) 1073b1a0a8a3SJed Brown { 10748f3b4b4dSDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 1075b1a0a8a3SJed Brown PetscErrorCode ierr; 1076b1a0a8a3SJed Brown 1077b1a0a8a3SJed Brown PetscFunctionBegin; 1078b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 10796a4f0f73SDmitry Karpeev ierr = PetscTryMethod(pc,"PCGASMSetSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,n,iis,ois));CHKERRQ(ierr); 10808f3b4b4dSDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 1081b1a0a8a3SJed Brown PetscFunctionReturn(0); 1082b1a0a8a3SJed Brown } 1083b1a0a8a3SJed Brown 1084b1a0a8a3SJed Brown 1085b1a0a8a3SJed Brown /*@ 1086f746d493SDmitry Karpeev PCGASMSetOverlap - Sets the overlap between a pair of subdomains for the 1087b1a0a8a3SJed Brown additive Schwarz preconditioner. Either all or no processors in the 10888f3b4b4dSDmitry Karpeev pc communicator must call this routine. 1089b1a0a8a3SJed Brown 10908f3b4b4dSDmitry Karpeev Logically Collective on pc 1091b1a0a8a3SJed Brown 1092b1a0a8a3SJed Brown Input Parameters: 1093b1a0a8a3SJed Brown + pc - the preconditioner context 10948f3b4b4dSDmitry Karpeev - ovl - the amount of overlap between subdomains (ovl >= 0, default value = 0) 1095b1a0a8a3SJed Brown 1096b1a0a8a3SJed Brown Options Database Key: 109706b43e7eSDmitry Karpeev . -pc_gasm_overlap <overlap> - Sets overlap 1098b1a0a8a3SJed Brown 1099b1a0a8a3SJed Brown Notes: 110006b43e7eSDmitry Karpeev By default the GASM preconditioner uses 1 subdomain per processor. To use 11018f3b4b4dSDmitry Karpeev multiple subdomain per perocessor or "straddling" subdomains that intersect 11028f3b4b4dSDmitry Karpeev multiple processors use PCGASMSetSubdomains() (or option -pc_gasm_total_subdomains <n>). 1103b1a0a8a3SJed Brown 11048f3b4b4dSDmitry Karpeev The overlap defaults to 0, so if one desires that no additional 1105b1a0a8a3SJed Brown overlap be computed beyond what may have been set with a call to 11068f3b4b4dSDmitry Karpeev PCGASMSetSubdomains(), then ovl must be set to be 0. In particular, if one does 11078f3b4b4dSDmitry Karpeev not explicitly set the subdomains in application code, then all overlap would be computed 1108f746d493SDmitry Karpeev internally by PETSc, and using an overlap of 0 would result in an GASM 1109b1a0a8a3SJed Brown variant that is equivalent to the block Jacobi preconditioner. 1110b1a0a8a3SJed Brown 1111b1a0a8a3SJed Brown Note that one can define initial index sets with any overlap via 111206b43e7eSDmitry Karpeev PCGASMSetSubdomains(); the routine PCGASMSetOverlap() merely allows 111306b43e7eSDmitry Karpeev PETSc to extend that overlap further, if desired. 1114b1a0a8a3SJed Brown 1115b1a0a8a3SJed Brown Level: intermediate 1116b1a0a8a3SJed Brown 1117f746d493SDmitry Karpeev .keywords: PC, GASM, set, overlap 1118b1a0a8a3SJed Brown 11198f3b4b4dSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMGetSubKSP(), 112006b43e7eSDmitry Karpeev PCGASMCreateSubdomains2D(), PCGASMGetSubdomains() 1121b1a0a8a3SJed Brown @*/ 11227087cfbeSBarry Smith PetscErrorCode PCGASMSetOverlap(PC pc,PetscInt ovl) 1123b1a0a8a3SJed Brown { 1124b1a0a8a3SJed Brown PetscErrorCode ierr; 11258f3b4b4dSDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 1126b1a0a8a3SJed Brown 1127b1a0a8a3SJed Brown PetscFunctionBegin; 1128b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1129b1a0a8a3SJed Brown PetscValidLogicalCollectiveInt(pc,ovl,2); 1130f746d493SDmitry Karpeev ierr = PetscTryMethod(pc,"PCGASMSetOverlap_C",(PC,PetscInt),(pc,ovl));CHKERRQ(ierr); 11318f3b4b4dSDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 1132b1a0a8a3SJed Brown PetscFunctionReturn(0); 1133b1a0a8a3SJed Brown } 1134b1a0a8a3SJed Brown 1135b1a0a8a3SJed Brown /*@ 1136f746d493SDmitry Karpeev PCGASMSetType - Sets the type of restriction and interpolation used 1137b1a0a8a3SJed Brown for local problems in the additive Schwarz method. 1138b1a0a8a3SJed Brown 1139b1a0a8a3SJed Brown Logically Collective on PC 1140b1a0a8a3SJed Brown 1141b1a0a8a3SJed Brown Input Parameters: 1142b1a0a8a3SJed Brown + pc - the preconditioner context 1143f746d493SDmitry Karpeev - type - variant of GASM, one of 1144b1a0a8a3SJed Brown .vb 1145f746d493SDmitry Karpeev PC_GASM_BASIC - full interpolation and restriction 1146f746d493SDmitry Karpeev PC_GASM_RESTRICT - full restriction, local processor interpolation 1147f746d493SDmitry Karpeev PC_GASM_INTERPOLATE - full interpolation, local processor restriction 1148f746d493SDmitry Karpeev PC_GASM_NONE - local processor restriction and interpolation 1149b1a0a8a3SJed Brown .ve 1150b1a0a8a3SJed Brown 1151b1a0a8a3SJed Brown Options Database Key: 1152f746d493SDmitry Karpeev . -pc_gasm_type [basic,restrict,interpolate,none] - Sets GASM type 1153b1a0a8a3SJed Brown 1154b1a0a8a3SJed Brown Level: intermediate 1155b1a0a8a3SJed Brown 1156f746d493SDmitry Karpeev .keywords: PC, GASM, set, type 1157b1a0a8a3SJed Brown 11588f3b4b4dSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMGetSubKSP(), 1159f746d493SDmitry Karpeev PCGASMCreateSubdomains2D() 1160b1a0a8a3SJed Brown @*/ 11617087cfbeSBarry Smith PetscErrorCode PCGASMSetType(PC pc,PCGASMType type) 1162b1a0a8a3SJed Brown { 1163b1a0a8a3SJed Brown PetscErrorCode ierr; 1164b1a0a8a3SJed Brown 1165b1a0a8a3SJed Brown PetscFunctionBegin; 1166b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1167b1a0a8a3SJed Brown PetscValidLogicalCollectiveEnum(pc,type,2); 1168f746d493SDmitry Karpeev ierr = PetscTryMethod(pc,"PCGASMSetType_C",(PC,PCGASMType),(pc,type));CHKERRQ(ierr); 1169b1a0a8a3SJed Brown PetscFunctionReturn(0); 1170b1a0a8a3SJed Brown } 1171b1a0a8a3SJed Brown 1172b1a0a8a3SJed Brown /*@ 1173f746d493SDmitry Karpeev PCGASMSetSortIndices - Determines whether subdomain indices are sorted. 1174b1a0a8a3SJed Brown 1175b1a0a8a3SJed Brown Logically Collective on PC 1176b1a0a8a3SJed Brown 1177b1a0a8a3SJed Brown Input Parameters: 1178b1a0a8a3SJed Brown + pc - the preconditioner context 1179b1a0a8a3SJed Brown - doSort - sort the subdomain indices 1180b1a0a8a3SJed Brown 1181b1a0a8a3SJed Brown Level: intermediate 1182b1a0a8a3SJed Brown 1183f746d493SDmitry Karpeev .keywords: PC, GASM, set, type 1184b1a0a8a3SJed Brown 11858f3b4b4dSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMGetSubKSP(), 1186f746d493SDmitry Karpeev PCGASMCreateSubdomains2D() 1187b1a0a8a3SJed Brown @*/ 11887087cfbeSBarry Smith PetscErrorCode PCGASMSetSortIndices(PC pc,PetscBool doSort) 1189b1a0a8a3SJed Brown { 1190b1a0a8a3SJed Brown PetscErrorCode ierr; 1191b1a0a8a3SJed Brown 1192b1a0a8a3SJed Brown PetscFunctionBegin; 1193b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1194b1a0a8a3SJed Brown PetscValidLogicalCollectiveBool(pc,doSort,2); 1195f746d493SDmitry Karpeev ierr = PetscTryMethod(pc,"PCGASMSetSortIndices_C",(PC,PetscBool),(pc,doSort));CHKERRQ(ierr); 1196b1a0a8a3SJed Brown PetscFunctionReturn(0); 1197b1a0a8a3SJed Brown } 1198b1a0a8a3SJed Brown 1199b1a0a8a3SJed Brown /*@C 1200f746d493SDmitry Karpeev PCGASMGetSubKSP - Gets the local KSP contexts for all blocks on 1201b1a0a8a3SJed Brown this processor. 1202b1a0a8a3SJed Brown 1203b1a0a8a3SJed Brown Collective on PC iff first_local is requested 1204b1a0a8a3SJed Brown 1205b1a0a8a3SJed Brown Input Parameter: 1206b1a0a8a3SJed Brown . pc - the preconditioner context 1207b1a0a8a3SJed Brown 1208b1a0a8a3SJed Brown Output Parameters: 12090298fd71SBarry Smith + n_local - the number of blocks on this processor or NULL 12100298fd71SBarry Smith . first_local - the global number of the first block on this processor or NULL, 12110298fd71SBarry Smith all processors must request or all must pass NULL 1212b1a0a8a3SJed Brown - ksp - the array of KSP contexts 1213b1a0a8a3SJed Brown 1214b1a0a8a3SJed Brown Note: 1215f746d493SDmitry Karpeev After PCGASMGetSubKSP() the array of KSPes is not to be freed 1216b1a0a8a3SJed Brown 1217b1a0a8a3SJed Brown Currently for some matrix implementations only 1 block per processor 1218b1a0a8a3SJed Brown is supported. 1219b1a0a8a3SJed Brown 1220f746d493SDmitry Karpeev You must call KSPSetUp() before calling PCGASMGetSubKSP(). 1221b1a0a8a3SJed Brown 1222b1a0a8a3SJed Brown Level: advanced 1223b1a0a8a3SJed Brown 1224f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, get, sub, KSP, context 1225b1a0a8a3SJed Brown 12268f3b4b4dSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMSetOverlap(), 1227f746d493SDmitry Karpeev PCGASMCreateSubdomains2D(), 1228b1a0a8a3SJed Brown @*/ 12297087cfbeSBarry Smith PetscErrorCode PCGASMGetSubKSP(PC pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp[]) 1230b1a0a8a3SJed Brown { 1231b1a0a8a3SJed Brown PetscErrorCode ierr; 1232b1a0a8a3SJed Brown 1233b1a0a8a3SJed Brown PetscFunctionBegin; 1234b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1235f746d493SDmitry Karpeev ierr = PetscUseMethod(pc,"PCGASMGetSubKSP_C",(PC,PetscInt*,PetscInt*,KSP **),(pc,n_local,first_local,ksp));CHKERRQ(ierr); 1236b1a0a8a3SJed Brown PetscFunctionReturn(0); 1237b1a0a8a3SJed Brown } 1238b1a0a8a3SJed Brown 1239b1a0a8a3SJed Brown /* -------------------------------------------------------------------------------------*/ 1240b1a0a8a3SJed Brown /*MC 1241f746d493SDmitry Karpeev PCGASM - Use the (restricted) additive Schwarz method, each block is (approximately) solved with 1242b1a0a8a3SJed Brown its own KSP object. 1243b1a0a8a3SJed Brown 1244b1a0a8a3SJed Brown Options Database Keys: 12458f3b4b4dSDmitry Karpeev + -pc_gasm_total_subdomains <n> - Sets total number of local subdomains to be distributed among processors 124606b43e7eSDmitry Karpeev . -pc_gasm_view_subdomains - activates the printing of subdomain indices in PCView(), -ksp_view or -snes_view 124706b43e7eSDmitry Karpeev . -pc_gasm_print_subdomains - activates the printing of subdomain indices in PCSetUp() 124806b43e7eSDmitry Karpeev . -pc_gasm_overlap <ovl> - Sets overlap by which to (automatically) extend local subdomains 1249f746d493SDmitry Karpeev - -pc_gasm_type [basic,restrict,interpolate,none] - Sets GASM type 1250b1a0a8a3SJed Brown 1251b1a0a8a3SJed Brown IMPORTANT: If you run with, for example, 3 blocks on 1 processor or 3 blocks on 3 processors you 1252f746d493SDmitry Karpeev will get a different convergence rate due to the default option of -pc_gasm_type restrict. Use 1253f746d493SDmitry Karpeev -pc_gasm_type basic to use the standard GASM. 1254b1a0a8a3SJed Brown 125595452b02SPatrick Sanan Notes: 125695452b02SPatrick Sanan Blocks can be shared by multiple processes. 1257b1a0a8a3SJed Brown 1258b1a0a8a3SJed Brown To set options on the solvers for each block append -sub_ to all the KSP, and PC 1259b1a0a8a3SJed Brown options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly 1260b1a0a8a3SJed Brown 1261f746d493SDmitry Karpeev To set the options on the solvers separate for each block call PCGASMGetSubKSP() 1262b1a0a8a3SJed Brown and set the options directly on the resulting KSP object (you can access its PC 1263b1a0a8a3SJed Brown with KSPGetPC()) 1264b1a0a8a3SJed Brown 1265b1a0a8a3SJed Brown 1266b1a0a8a3SJed Brown Level: beginner 1267b1a0a8a3SJed Brown 1268b1a0a8a3SJed Brown Concepts: additive Schwarz method 1269b1a0a8a3SJed Brown 1270b1a0a8a3SJed Brown References: 127196a0c994SBarry Smith + 1. - M Dryja, OB Widlund, An additive variant of the Schwarz alternating method for the case of many subregions 127296a0c994SBarry Smith Courant Institute, New York University Technical report 127396a0c994SBarry Smith - 2. - Barry Smith, Petter Bjorstad, and William Gropp, Domain Decompositions: Parallel Multilevel Methods for Elliptic Partial Differential Equations, 127496a0c994SBarry Smith Cambridge University Press. 1275b1a0a8a3SJed Brown 1276b1a0a8a3SJed Brown .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 127749517cdeSBarry Smith PCBJACOBI, PCGASMGetSubKSP(), PCGASMSetSubdomains(), 12788f3b4b4dSDmitry Karpeev PCSetModifySubmatrices(), PCGASMSetOverlap(), PCGASMSetType() 1279b1a0a8a3SJed Brown 1280b1a0a8a3SJed Brown M*/ 1281b1a0a8a3SJed Brown 12828cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_GASM(PC pc) 1283b1a0a8a3SJed Brown { 1284b1a0a8a3SJed Brown PetscErrorCode ierr; 1285f746d493SDmitry Karpeev PC_GASM *osm; 1286b1a0a8a3SJed Brown 1287b1a0a8a3SJed Brown PetscFunctionBegin; 1288b00a9115SJed Brown ierr = PetscNewLog(pc,&osm);CHKERRQ(ierr); 12892fa5cd67SKarl Rupp 12908f3b4b4dSDmitry Karpeev osm->N = PETSC_DETERMINE; 129106b43e7eSDmitry Karpeev osm->n = PETSC_DECIDE; 12928f3b4b4dSDmitry Karpeev osm->nmax = PETSC_DETERMINE; 12938f3b4b4dSDmitry Karpeev osm->overlap = 0; 1294b1a0a8a3SJed Brown osm->ksp = 0; 12956a4f0f73SDmitry Karpeev osm->gorestriction = 0; 12966a4f0f73SDmitry Karpeev osm->girestriction = 0; 1297ea91fabdSFande Kong osm->pctoouter = 0; 1298ab3e923fSDmitry Karpeev osm->gx = 0; 1299ab3e923fSDmitry Karpeev osm->gy = 0; 1300b1a0a8a3SJed Brown osm->x = 0; 1301b1a0a8a3SJed Brown osm->y = 0; 1302ea91fabdSFande Kong osm->pcx = 0; 1303ea91fabdSFande Kong osm->pcy = 0; 1304ea91fabdSFande Kong osm->permutationIS = 0; 1305ea91fabdSFande Kong osm->permutationP = 0; 1306ea91fabdSFande Kong osm->pcmat = 0; 13076a4f0f73SDmitry Karpeev osm->ois = 0; 13086a4f0f73SDmitry Karpeev osm->iis = 0; 1309b1a0a8a3SJed Brown osm->pmat = 0; 1310f746d493SDmitry Karpeev osm->type = PC_GASM_RESTRICT; 13116a4f0f73SDmitry Karpeev osm->same_subdomain_solvers = PETSC_TRUE; 1312b1a0a8a3SJed Brown osm->sort_indices = PETSC_TRUE; 1313d709ea83SDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 1314ea91fabdSFande Kong osm->hierarchicalpartitioning = PETSC_FALSE; 1315b1a0a8a3SJed Brown 1316b1a0a8a3SJed Brown pc->data = (void*)osm; 1317f746d493SDmitry Karpeev pc->ops->apply = PCApply_GASM; 1318f746d493SDmitry Karpeev pc->ops->applytranspose = PCApplyTranspose_GASM; 1319f746d493SDmitry Karpeev pc->ops->setup = PCSetUp_GASM; 1320a06653b4SBarry Smith pc->ops->reset = PCReset_GASM; 1321f746d493SDmitry Karpeev pc->ops->destroy = PCDestroy_GASM; 1322f746d493SDmitry Karpeev pc->ops->setfromoptions = PCSetFromOptions_GASM; 1323f746d493SDmitry Karpeev pc->ops->setuponblocks = PCSetUpOnBlocks_GASM; 1324f746d493SDmitry Karpeev pc->ops->view = PCView_GASM; 1325b1a0a8a3SJed Brown pc->ops->applyrichardson = 0; 1326b1a0a8a3SJed Brown 1327bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetSubdomains_C",PCGASMSetSubdomains_GASM);CHKERRQ(ierr); 1328bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetOverlap_C",PCGASMSetOverlap_GASM);CHKERRQ(ierr); 1329bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetType_C",PCGASMSetType_GASM);CHKERRQ(ierr); 1330bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetSortIndices_C",PCGASMSetSortIndices_GASM);CHKERRQ(ierr); 1331bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCGASMGetSubKSP_C",PCGASMGetSubKSP_GASM);CHKERRQ(ierr); 1332b1a0a8a3SJed Brown PetscFunctionReturn(0); 1333b1a0a8a3SJed Brown } 1334b1a0a8a3SJed Brown 1335b1a0a8a3SJed Brown 13368f3b4b4dSDmitry Karpeev PetscErrorCode PCGASMCreateLocalSubdomains(Mat A, PetscInt nloc, IS *iis[]) 1337b1a0a8a3SJed Brown { 1338b1a0a8a3SJed Brown MatPartitioning mpart; 1339b1a0a8a3SJed Brown const char *prefix; 1340b1a0a8a3SJed Brown PetscInt i,j,rstart,rend,bs; 1341976e8c5aSLisandro Dalcin PetscBool hasop, isbaij = PETSC_FALSE,foundpart = PETSC_FALSE; 13420298fd71SBarry Smith Mat Ad = NULL, adj; 1343b1a0a8a3SJed Brown IS ispart,isnumb,*is; 1344b1a0a8a3SJed Brown PetscErrorCode ierr; 1345b1a0a8a3SJed Brown 1346b1a0a8a3SJed Brown PetscFunctionBegin; 13478f3b4b4dSDmitry Karpeev if (nloc < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of local subdomains must > 0, got nloc = %D",nloc); 1348b1a0a8a3SJed Brown 1349b1a0a8a3SJed Brown /* Get prefix, row distribution, and block size */ 1350b1a0a8a3SJed Brown ierr = MatGetOptionsPrefix(A,&prefix);CHKERRQ(ierr); 1351b1a0a8a3SJed Brown ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 1352b1a0a8a3SJed Brown ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr); 1353b1a0a8a3SJed 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); 1354b1a0a8a3SJed Brown 1355b1a0a8a3SJed Brown /* Get diagonal block from matrix if possible */ 1356976e8c5aSLisandro Dalcin ierr = MatHasOperation(A,MATOP_GET_DIAGONAL_BLOCK,&hasop);CHKERRQ(ierr); 1357976e8c5aSLisandro Dalcin if (hasop) { 135811bd1e4dSLisandro Dalcin ierr = MatGetDiagonalBlock(A,&Ad);CHKERRQ(ierr); 1359b1a0a8a3SJed Brown } 1360b1a0a8a3SJed Brown if (Ad) { 1361251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)Ad,MATSEQBAIJ,&isbaij);CHKERRQ(ierr); 1362251f4c67SDmitry Karpeev if (!isbaij) {ierr = PetscObjectTypeCompare((PetscObject)Ad,MATSEQSBAIJ,&isbaij);CHKERRQ(ierr);} 1363b1a0a8a3SJed Brown } 13648f3b4b4dSDmitry Karpeev if (Ad && nloc > 1) { 1365b1a0a8a3SJed Brown PetscBool match,done; 1366b1a0a8a3SJed Brown /* Try to setup a good matrix partitioning if available */ 1367b1a0a8a3SJed Brown ierr = MatPartitioningCreate(PETSC_COMM_SELF,&mpart);CHKERRQ(ierr); 1368b1a0a8a3SJed Brown ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr); 1369b1a0a8a3SJed Brown ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr); 1370251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGCURRENT,&match);CHKERRQ(ierr); 1371b1a0a8a3SJed Brown if (!match) { 1372251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGSQUARE,&match);CHKERRQ(ierr); 1373b1a0a8a3SJed Brown } 1374b1a0a8a3SJed Brown if (!match) { /* assume a "good" partitioner is available */ 13751a83f524SJed Brown PetscInt na; 13761a83f524SJed Brown const PetscInt *ia,*ja; 1377b1a0a8a3SJed Brown ierr = MatGetRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr); 1378b1a0a8a3SJed Brown if (done) { 1379b1a0a8a3SJed Brown /* Build adjacency matrix by hand. Unfortunately a call to 1380b1a0a8a3SJed Brown MatConvert(Ad,MATMPIADJ,MAT_INITIAL_MATRIX,&adj) will 1381b1a0a8a3SJed Brown remove the block-aij structure and we cannot expect 1382b1a0a8a3SJed Brown MatPartitioning to split vertices as we need */ 13831a83f524SJed Brown PetscInt i,j,len,nnz,cnt,*iia=0,*jja=0; 13841a83f524SJed Brown const PetscInt *row; 1385b1a0a8a3SJed Brown nnz = 0; 1386b1a0a8a3SJed Brown for (i=0; i<na; i++) { /* count number of nonzeros */ 1387b1a0a8a3SJed Brown len = ia[i+1] - ia[i]; 1388b1a0a8a3SJed Brown row = ja + ia[i]; 1389b1a0a8a3SJed Brown for (j=0; j<len; j++) { 1390b1a0a8a3SJed Brown if (row[j] == i) { /* don't count diagonal */ 1391b1a0a8a3SJed Brown len--; break; 1392b1a0a8a3SJed Brown } 1393b1a0a8a3SJed Brown } 1394b1a0a8a3SJed Brown nnz += len; 1395b1a0a8a3SJed Brown } 1396854ce69bSBarry Smith ierr = PetscMalloc1(na+1,&iia);CHKERRQ(ierr); 1397854ce69bSBarry Smith ierr = PetscMalloc1(nnz,&jja);CHKERRQ(ierr); 1398b1a0a8a3SJed Brown nnz = 0; 1399b1a0a8a3SJed Brown iia[0] = 0; 1400b1a0a8a3SJed Brown for (i=0; i<na; i++) { /* fill adjacency */ 1401b1a0a8a3SJed Brown cnt = 0; 1402b1a0a8a3SJed Brown len = ia[i+1] - ia[i]; 1403b1a0a8a3SJed Brown row = ja + ia[i]; 1404b1a0a8a3SJed Brown for (j=0; j<len; j++) { 14052fa5cd67SKarl Rupp if (row[j] != i) jja[nnz+cnt++] = row[j]; /* if not diagonal */ 1406b1a0a8a3SJed Brown } 1407b1a0a8a3SJed Brown nnz += cnt; 1408b1a0a8a3SJed Brown iia[i+1] = nnz; 1409b1a0a8a3SJed Brown } 1410b1a0a8a3SJed Brown /* Partitioning of the adjacency matrix */ 14110298fd71SBarry Smith ierr = MatCreateMPIAdj(PETSC_COMM_SELF,na,na,iia,jja,NULL,&adj);CHKERRQ(ierr); 1412b1a0a8a3SJed Brown ierr = MatPartitioningSetAdjacency(mpart,adj);CHKERRQ(ierr); 14138f3b4b4dSDmitry Karpeev ierr = MatPartitioningSetNParts(mpart,nloc);CHKERRQ(ierr); 1414b1a0a8a3SJed Brown ierr = MatPartitioningApply(mpart,&ispart);CHKERRQ(ierr); 1415b1a0a8a3SJed Brown ierr = ISPartitioningToNumbering(ispart,&isnumb);CHKERRQ(ierr); 14166bf464f9SBarry Smith ierr = MatDestroy(&adj);CHKERRQ(ierr); 1417b1a0a8a3SJed Brown foundpart = PETSC_TRUE; 1418b1a0a8a3SJed Brown } 1419b1a0a8a3SJed Brown ierr = MatRestoreRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr); 1420b1a0a8a3SJed Brown } 14216bf464f9SBarry Smith ierr = MatPartitioningDestroy(&mpart);CHKERRQ(ierr); 1422b1a0a8a3SJed Brown } 14238f3b4b4dSDmitry Karpeev ierr = PetscMalloc1(nloc,&is);CHKERRQ(ierr); 1424b1a0a8a3SJed Brown if (!foundpart) { 1425b1a0a8a3SJed Brown 1426b1a0a8a3SJed Brown /* Partitioning by contiguous chunks of rows */ 1427b1a0a8a3SJed Brown 1428b1a0a8a3SJed Brown PetscInt mbs = (rend-rstart)/bs; 1429b1a0a8a3SJed Brown PetscInt start = rstart; 14308f3b4b4dSDmitry Karpeev for (i=0; i<nloc; i++) { 14318f3b4b4dSDmitry Karpeev PetscInt count = (mbs/nloc + ((mbs % nloc) > i)) * bs; 1432b1a0a8a3SJed Brown ierr = ISCreateStride(PETSC_COMM_SELF,count,start,1,&is[i]);CHKERRQ(ierr); 1433b1a0a8a3SJed Brown start += count; 1434b1a0a8a3SJed Brown } 1435b1a0a8a3SJed Brown 1436b1a0a8a3SJed Brown } else { 1437b1a0a8a3SJed Brown 1438b1a0a8a3SJed Brown /* Partitioning by adjacency of diagonal block */ 1439b1a0a8a3SJed Brown 1440b1a0a8a3SJed Brown const PetscInt *numbering; 1441b1a0a8a3SJed Brown PetscInt *count,nidx,*indices,*newidx,start=0; 1442b1a0a8a3SJed Brown /* Get node count in each partition */ 14438f3b4b4dSDmitry Karpeev ierr = PetscMalloc1(nloc,&count);CHKERRQ(ierr); 14448f3b4b4dSDmitry Karpeev ierr = ISPartitioningCount(ispart,nloc,count);CHKERRQ(ierr); 1445b1a0a8a3SJed Brown if (isbaij && bs > 1) { /* adjust for the block-aij case */ 14468f3b4b4dSDmitry Karpeev for (i=0; i<nloc; i++) count[i] *= bs; 1447b1a0a8a3SJed Brown } 1448b1a0a8a3SJed Brown /* Build indices from node numbering */ 1449b1a0a8a3SJed Brown ierr = ISGetLocalSize(isnumb,&nidx);CHKERRQ(ierr); 1450785e854fSJed Brown ierr = PetscMalloc1(nidx,&indices);CHKERRQ(ierr); 1451b1a0a8a3SJed Brown for (i=0; i<nidx; i++) indices[i] = i; /* needs to be initialized */ 1452b1a0a8a3SJed Brown ierr = ISGetIndices(isnumb,&numbering);CHKERRQ(ierr); 1453b1a0a8a3SJed Brown ierr = PetscSortIntWithPermutation(nidx,numbering,indices);CHKERRQ(ierr); 1454b1a0a8a3SJed Brown ierr = ISRestoreIndices(isnumb,&numbering);CHKERRQ(ierr); 1455b1a0a8a3SJed Brown if (isbaij && bs > 1) { /* adjust for the block-aij case */ 1456785e854fSJed Brown ierr = PetscMalloc1(nidx*bs,&newidx);CHKERRQ(ierr); 14572fa5cd67SKarl Rupp for (i=0; i<nidx; i++) { 14582fa5cd67SKarl Rupp for (j=0; j<bs; j++) newidx[i*bs+j] = indices[i]*bs + j; 14592fa5cd67SKarl Rupp } 1460b1a0a8a3SJed Brown ierr = PetscFree(indices);CHKERRQ(ierr); 1461b1a0a8a3SJed Brown nidx *= bs; 1462b1a0a8a3SJed Brown indices = newidx; 1463b1a0a8a3SJed Brown } 1464b1a0a8a3SJed Brown /* Shift to get global indices */ 1465b1a0a8a3SJed Brown for (i=0; i<nidx; i++) indices[i] += rstart; 1466b1a0a8a3SJed Brown 1467b1a0a8a3SJed Brown /* Build the index sets for each block */ 14688f3b4b4dSDmitry Karpeev for (i=0; i<nloc; i++) { 1469b1a0a8a3SJed Brown ierr = ISCreateGeneral(PETSC_COMM_SELF,count[i],&indices[start],PETSC_COPY_VALUES,&is[i]);CHKERRQ(ierr); 1470b1a0a8a3SJed Brown ierr = ISSort(is[i]);CHKERRQ(ierr); 1471b1a0a8a3SJed Brown start += count[i]; 1472b1a0a8a3SJed Brown } 1473b1a0a8a3SJed Brown 1474302440fdSBarry Smith ierr = PetscFree(count);CHKERRQ(ierr); 1475302440fdSBarry Smith ierr = PetscFree(indices);CHKERRQ(ierr); 1476fcfd50ebSBarry Smith ierr = ISDestroy(&isnumb);CHKERRQ(ierr); 1477fcfd50ebSBarry Smith ierr = ISDestroy(&ispart);CHKERRQ(ierr); 1478b1a0a8a3SJed Brown } 14796a4f0f73SDmitry Karpeev *iis = is; 14808f3b4b4dSDmitry Karpeev PetscFunctionReturn(0); 14818f3b4b4dSDmitry Karpeev } 14828f3b4b4dSDmitry Karpeev 1483b541e6a4SDmitry Karpeev PETSC_INTERN PetscErrorCode PCGASMCreateStraddlingSubdomains(Mat A,PetscInt N,PetscInt *n,IS *iis[]) 14848f3b4b4dSDmitry Karpeev { 14858f3b4b4dSDmitry Karpeev PetscErrorCode ierr; 14868f3b4b4dSDmitry Karpeev 14878f3b4b4dSDmitry Karpeev PetscFunctionBegin; 14888f3b4b4dSDmitry Karpeev ierr = MatSubdomainsCreateCoalesce(A,N,n,iis);CHKERRQ(ierr); 14898f3b4b4dSDmitry Karpeev PetscFunctionReturn(0); 14908f3b4b4dSDmitry Karpeev } 14918f3b4b4dSDmitry Karpeev 14928f3b4b4dSDmitry Karpeev 14938f3b4b4dSDmitry Karpeev 14948f3b4b4dSDmitry Karpeev /*@C 14958f3b4b4dSDmitry Karpeev PCGASMCreateSubdomains - Creates n index sets defining n nonoverlapping subdomains for the additive 14968f3b4b4dSDmitry Karpeev Schwarz preconditioner for a any problem based on its matrix. 14978f3b4b4dSDmitry Karpeev 14988f3b4b4dSDmitry Karpeev Collective 14998f3b4b4dSDmitry Karpeev 15008f3b4b4dSDmitry Karpeev Input Parameters: 15018f3b4b4dSDmitry Karpeev + A - The global matrix operator 15028f3b4b4dSDmitry Karpeev - N - the number of global subdomains requested 15038f3b4b4dSDmitry Karpeev 15048f3b4b4dSDmitry Karpeev Output Parameters: 15058f3b4b4dSDmitry Karpeev + n - the number of subdomains created on this processor 15068f3b4b4dSDmitry Karpeev - iis - the array of index sets defining the local inner subdomains (on which the correction is applied) 15078f3b4b4dSDmitry Karpeev 15088f3b4b4dSDmitry Karpeev Level: advanced 15098f3b4b4dSDmitry Karpeev 15108f3b4b4dSDmitry Karpeev Note: When N >= A's communicator size, each subdomain is local -- contained within a single processor. 15118f3b4b4dSDmitry Karpeev When N < size, the subdomains are 'straddling' (processor boundaries) and are no longer local. 15128f3b4b4dSDmitry Karpeev The resulting subdomains can be use in PCGASMSetSubdomains(pc,n,iss,NULL). The overlapping 15138f3b4b4dSDmitry Karpeev outer subdomains will be automatically generated from these according to the requested amount of 15148f3b4b4dSDmitry Karpeev overlap; this is currently supported only with local subdomains. 15158f3b4b4dSDmitry Karpeev 15168f3b4b4dSDmitry Karpeev 15178f3b4b4dSDmitry Karpeev .keywords: PC, GASM, additive Schwarz, create, subdomains, unstructured grid 15188f3b4b4dSDmitry Karpeev 15198f3b4b4dSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMDestroySubdomains() 15208f3b4b4dSDmitry Karpeev @*/ 1521b541e6a4SDmitry Karpeev PetscErrorCode PCGASMCreateSubdomains(Mat A,PetscInt N,PetscInt *n,IS *iis[]) 15228f3b4b4dSDmitry Karpeev { 15238f3b4b4dSDmitry Karpeev PetscMPIInt size; 15248f3b4b4dSDmitry Karpeev PetscErrorCode ierr; 15258f3b4b4dSDmitry Karpeev 15268f3b4b4dSDmitry Karpeev PetscFunctionBegin; 15278f3b4b4dSDmitry Karpeev PetscValidHeaderSpecific(A,MAT_CLASSID,1); 15288f3b4b4dSDmitry Karpeev PetscValidPointer(iis,4); 15298f3b4b4dSDmitry Karpeev 15308f3b4b4dSDmitry Karpeev if (N < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Number of subdomains must be > 0, N = %D",N); 15318f3b4b4dSDmitry Karpeev ierr = MPI_Comm_size(PetscObjectComm((PetscObject)A),&size);CHKERRQ(ierr); 15328f3b4b4dSDmitry Karpeev if (N >= size) { 15338f3b4b4dSDmitry Karpeev *n = N/size + (N%size); 15348f3b4b4dSDmitry Karpeev ierr = PCGASMCreateLocalSubdomains(A,*n,iis);CHKERRQ(ierr); 15356a4f0f73SDmitry Karpeev } else { 15368f3b4b4dSDmitry Karpeev ierr = PCGASMCreateStraddlingSubdomains(A,N,n,iis);CHKERRQ(ierr); 15376a4f0f73SDmitry Karpeev } 1538b1a0a8a3SJed Brown PetscFunctionReturn(0); 1539b1a0a8a3SJed Brown } 1540b1a0a8a3SJed Brown 1541b1a0a8a3SJed Brown /*@C 1542f746d493SDmitry Karpeev PCGASMDestroySubdomains - Destroys the index sets created with 15438f3b4b4dSDmitry Karpeev PCGASMCreateSubdomains() or PCGASMCreateSubdomains2D. Should be 154406b43e7eSDmitry Karpeev called after setting subdomains with PCGASMSetSubdomains(). 1545b1a0a8a3SJed Brown 1546b1a0a8a3SJed Brown Collective 1547b1a0a8a3SJed Brown 1548b1a0a8a3SJed Brown Input Parameters: 1549b1a0a8a3SJed Brown + n - the number of index sets 15506a4f0f73SDmitry Karpeev . iis - the array of inner subdomains, 15510298fd71SBarry Smith - ois - the array of outer subdomains, can be NULL 1552b1a0a8a3SJed Brown 15536a4f0f73SDmitry Karpeev Level: intermediate 15546a4f0f73SDmitry Karpeev 155595452b02SPatrick Sanan Notes: 155695452b02SPatrick Sanan this is merely a convenience subroutine that walks each list, 15572c112581SDmitry Karpeev destroys each IS on the list, and then frees the list. At the end the 15582c112581SDmitry Karpeev list pointers are set to NULL. 1559b1a0a8a3SJed Brown 1560f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, create, subdomains, unstructured grid 1561b1a0a8a3SJed Brown 15628f3b4b4dSDmitry Karpeev .seealso: PCGASMCreateSubdomains(), PCGASMSetSubdomains() 1563b1a0a8a3SJed Brown @*/ 15642c112581SDmitry Karpeev PetscErrorCode PCGASMDestroySubdomains(PetscInt n,IS **iis,IS **ois) 1565b1a0a8a3SJed Brown { 1566b1a0a8a3SJed Brown PetscInt i; 1567b1a0a8a3SJed Brown PetscErrorCode ierr; 15685fd66863SKarl Rupp 1569b1a0a8a3SJed Brown PetscFunctionBegin; 1570a35b7b57SDmitry Karpeev if (n <= 0) PetscFunctionReturn(0); 15716a4f0f73SDmitry Karpeev if (ois) { 15722c112581SDmitry Karpeev PetscValidPointer(ois,3); 15732c112581SDmitry Karpeev if (*ois) { 15742c112581SDmitry Karpeev PetscValidPointer(*ois,3); 1575a35b7b57SDmitry Karpeev for (i=0; i<n; i++) { 15762c112581SDmitry Karpeev ierr = ISDestroy(&(*ois)[i]);CHKERRQ(ierr); 1577a35b7b57SDmitry Karpeev } 1578135757f6SDmitry Karpeev ierr = PetscFree((*ois));CHKERRQ(ierr); 15792c112581SDmitry Karpeev } 1580b1a0a8a3SJed Brown } 1581b9d0fdaaSFande Kong if (iis) { 1582b9d0fdaaSFande Kong PetscValidPointer(iis,2); 1583b9d0fdaaSFande Kong if (*iis) { 1584b9d0fdaaSFande Kong PetscValidPointer(*iis,2); 1585b9d0fdaaSFande Kong for (i=0; i<n; i++) { 1586b9d0fdaaSFande Kong ierr = ISDestroy(&(*iis)[i]);CHKERRQ(ierr); 1587b9d0fdaaSFande Kong } 1588b9d0fdaaSFande Kong ierr = PetscFree((*iis));CHKERRQ(ierr); 1589b9d0fdaaSFande Kong } 1590b9d0fdaaSFande Kong } 1591b1a0a8a3SJed Brown PetscFunctionReturn(0); 1592b1a0a8a3SJed Brown } 1593b1a0a8a3SJed Brown 1594af538404SDmitry Karpeev 1595af538404SDmitry Karpeev #define PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,xleft_loc,ylow_loc,xright_loc,yhigh_loc,n) \ 1596af538404SDmitry Karpeev { \ 1597af538404SDmitry Karpeev PetscInt first_row = first/M, last_row = last/M+1; \ 1598af538404SDmitry Karpeev /* \ 1599af538404SDmitry Karpeev Compute ylow_loc and yhigh_loc so that (ylow_loc,xleft) and (yhigh_loc,xright) are the corners \ 1600af538404SDmitry Karpeev of the bounding box of the intersection of the subdomain with the local ownership range (local \ 1601af538404SDmitry Karpeev subdomain). \ 1602af538404SDmitry Karpeev Also compute xleft_loc and xright_loc as the lower and upper bounds on the first and last rows \ 1603af538404SDmitry Karpeev of the intersection. \ 1604af538404SDmitry Karpeev */ \ 1605af538404SDmitry Karpeev /* ylow_loc is the grid row containing the first element of the local sumbdomain */ \ 1606eec7e3faSDmitry Karpeev *ylow_loc = PetscMax(first_row,ylow); \ 1607af538404SDmitry Karpeev /* xleft_loc is the offset of first element of the local subdomain within its grid row (might actually be outside the local subdomain) */ \ 1608eec7e3faSDmitry Karpeev *xleft_loc = *ylow_loc==first_row ? PetscMax(first%M,xleft) : xleft; \ 1609af538404SDmitry Karpeev /* yhigh_loc is the grid row above the last local subdomain element */ \ 1610eec7e3faSDmitry Karpeev *yhigh_loc = PetscMin(last_row,yhigh); \ 1611af538404SDmitry Karpeev /* xright is the offset of the end of the local subdomain within its grid row (might actually be outside the local subdomain) */ \ 1612eec7e3faSDmitry Karpeev *xright_loc = *yhigh_loc==last_row ? PetscMin(xright,last%M) : xright; \ 1613af538404SDmitry Karpeev /* Now compute the size of the local subdomain n. */ \ 1614c3518bceSDmitry Karpeev *n = 0; \ 1615eec7e3faSDmitry Karpeev if (*ylow_loc < *yhigh_loc) { \ 1616af538404SDmitry Karpeev PetscInt width = xright-xleft; \ 1617eec7e3faSDmitry Karpeev *n += width*(*yhigh_loc-*ylow_loc-1); \ 1618eec7e3faSDmitry Karpeev *n += PetscMin(PetscMax(*xright_loc-xleft,0),width); \ 1619eec7e3faSDmitry Karpeev *n -= PetscMin(PetscMax(*xleft_loc-xleft,0), width); \ 1620af538404SDmitry Karpeev } \ 1621af538404SDmitry Karpeev } 1622af538404SDmitry Karpeev 1623c3518bceSDmitry Karpeev 1624eec7e3faSDmitry Karpeev 1625b1a0a8a3SJed Brown /*@ 1626f746d493SDmitry Karpeev PCGASMCreateSubdomains2D - Creates the index sets for the overlapping Schwarz 1627b1a0a8a3SJed Brown preconditioner for a two-dimensional problem on a regular grid. 1628b1a0a8a3SJed Brown 1629af538404SDmitry Karpeev Collective 1630b1a0a8a3SJed Brown 1631b1a0a8a3SJed Brown Input Parameters: 163206b43e7eSDmitry Karpeev + M, N - the global number of grid points in the x and y directions 1633af538404SDmitry Karpeev . Mdomains, Ndomains - the global number of subdomains in the x and y directions 1634b1a0a8a3SJed Brown . dof - degrees of freedom per node 1635b1a0a8a3SJed Brown - overlap - overlap in mesh lines 1636b1a0a8a3SJed Brown 1637b1a0a8a3SJed Brown Output Parameters: 1638af538404SDmitry Karpeev + Nsub - the number of local subdomains created 16396a4f0f73SDmitry Karpeev . iis - array of index sets defining inner (nonoverlapping) subdomains 16406a4f0f73SDmitry Karpeev - ois - array of index sets defining outer (overlapping, if overlap > 0) subdomains 1641b1a0a8a3SJed Brown 1642b1a0a8a3SJed Brown 1643b1a0a8a3SJed Brown Level: advanced 1644b1a0a8a3SJed Brown 1645f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, create, subdomains, 2D, regular grid 1646b1a0a8a3SJed Brown 16478f3b4b4dSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMGetSubKSP(), PCGASMSetOverlap() 1648b1a0a8a3SJed Brown @*/ 16496a4f0f73SDmitry Karpeev PetscErrorCode PCGASMCreateSubdomains2D(PC pc,PetscInt M,PetscInt N,PetscInt Mdomains,PetscInt Ndomains,PetscInt dof,PetscInt overlap,PetscInt *nsub,IS **iis,IS **ois) 1650b1a0a8a3SJed Brown { 1651b1a0a8a3SJed Brown PetscErrorCode ierr; 16522bbb417fSDmitry Karpeev PetscMPIInt size, rank; 16532bbb417fSDmitry Karpeev PetscInt i, j; 16542bbb417fSDmitry Karpeev PetscInt maxheight, maxwidth; 16552bbb417fSDmitry Karpeev PetscInt xstart, xleft, xright, xleft_loc, xright_loc; 16562bbb417fSDmitry Karpeev PetscInt ystart, ylow, yhigh, ylow_loc, yhigh_loc; 1657eec7e3faSDmitry Karpeev PetscInt x[2][2], y[2][2], n[2]; 16582bbb417fSDmitry Karpeev PetscInt first, last; 16592bbb417fSDmitry Karpeev PetscInt nidx, *idx; 16602bbb417fSDmitry Karpeev PetscInt ii,jj,s,q,d; 1661af538404SDmitry Karpeev PetscInt k,kk; 16622bbb417fSDmitry Karpeev PetscMPIInt color; 16632bbb417fSDmitry Karpeev MPI_Comm comm, subcomm; 16646a4f0f73SDmitry Karpeev IS **xis = 0, **is = ois, **is_local = iis; 1665d34fcf5fSBarry Smith 1666b1a0a8a3SJed Brown PetscFunctionBegin; 16672bbb417fSDmitry Karpeev ierr = PetscObjectGetComm((PetscObject)pc, &comm);CHKERRQ(ierr); 16682bbb417fSDmitry Karpeev ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 16692bbb417fSDmitry Karpeev ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 16702bbb417fSDmitry Karpeev ierr = MatGetOwnershipRange(pc->pmat, &first, &last);CHKERRQ(ierr); 1671d34fcf5fSBarry 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) " 16722bbb417fSDmitry Karpeev "does not respect the number of degrees of freedom per grid point %D", first, last, dof); 1673d34fcf5fSBarry Smith 1674af538404SDmitry Karpeev /* Determine the number of domains with nonzero intersections with the local ownership range. */ 16752bbb417fSDmitry Karpeev s = 0; 1676af538404SDmitry Karpeev ystart = 0; 1677af538404SDmitry Karpeev for (j=0; j<Ndomains; ++j) { 1678af538404SDmitry Karpeev maxheight = N/Ndomains + ((N % Ndomains) > j); /* Maximal height of subdomain */ 1679af538404SDmitry Karpeev if (maxheight < 2) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many %D subdomains in the vertical directon for mesh height %D", Ndomains, N); 1680eec7e3faSDmitry Karpeev /* Vertical domain limits with an overlap. */ 1681eec7e3faSDmitry Karpeev ylow = PetscMax(ystart - overlap,0); 1682eec7e3faSDmitry Karpeev yhigh = PetscMin(ystart + maxheight + overlap,N); 1683b1a0a8a3SJed Brown xstart = 0; 1684af538404SDmitry Karpeev for (i=0; i<Mdomains; ++i) { 1685af538404SDmitry Karpeev maxwidth = M/Mdomains + ((M % Mdomains) > i); /* Maximal width of subdomain */ 1686af538404SDmitry 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); 1687eec7e3faSDmitry Karpeev /* Horizontal domain limits with an overlap. */ 1688eec7e3faSDmitry Karpeev xleft = PetscMax(xstart - overlap,0); 1689eec7e3faSDmitry Karpeev xright = PetscMin(xstart + maxwidth + overlap,M); 1690af538404SDmitry Karpeev /* 1691af538404SDmitry Karpeev Determine whether this subdomain intersects this processor's ownership range of pc->pmat. 1692af538404SDmitry Karpeev */ 1693c3518bceSDmitry Karpeev PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,(&xleft_loc),(&ylow_loc),(&xright_loc),(&yhigh_loc),(&nidx)); 16942fa5cd67SKarl Rupp if (nidx) ++s; 1695af538404SDmitry Karpeev xstart += maxwidth; 1696af538404SDmitry Karpeev } /* for (i = 0; i < Mdomains; ++i) */ 1697af538404SDmitry Karpeev ystart += maxheight; 1698af538404SDmitry Karpeev } /* for (j = 0; j < Ndomains; ++j) */ 16992fa5cd67SKarl Rupp 1700af538404SDmitry Karpeev /* Now we can allocate the necessary number of ISs. */ 1701af538404SDmitry Karpeev *nsub = s; 1702854ce69bSBarry Smith ierr = PetscMalloc1(*nsub,is);CHKERRQ(ierr); 1703854ce69bSBarry Smith ierr = PetscMalloc1(*nsub,is_local);CHKERRQ(ierr); 1704af538404SDmitry Karpeev s = 0; 1705af538404SDmitry Karpeev ystart = 0; 1706af538404SDmitry Karpeev for (j=0; j<Ndomains; ++j) { 1707af538404SDmitry Karpeev maxheight = N/Ndomains + ((N % Ndomains) > j); /* Maximal height of subdomain */ 1708af538404SDmitry 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); 1709af538404SDmitry Karpeev /* Vertical domain limits with an overlap. */ 1710eec7e3faSDmitry Karpeev y[0][0] = PetscMax(ystart - overlap,0); 1711eec7e3faSDmitry Karpeev y[0][1] = PetscMin(ystart + maxheight + overlap,N); 1712eec7e3faSDmitry Karpeev /* Vertical domain limits without an overlap. */ 1713eec7e3faSDmitry Karpeev y[1][0] = ystart; 1714eec7e3faSDmitry Karpeev y[1][1] = PetscMin(ystart + maxheight,N); 1715eec7e3faSDmitry Karpeev xstart = 0; 1716af538404SDmitry Karpeev for (i=0; i<Mdomains; ++i) { 1717af538404SDmitry Karpeev maxwidth = M/Mdomains + ((M % Mdomains) > i); /* Maximal width of subdomain */ 1718af538404SDmitry 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); 1719af538404SDmitry Karpeev /* Horizontal domain limits with an overlap. */ 1720eec7e3faSDmitry Karpeev x[0][0] = PetscMax(xstart - overlap,0); 1721eec7e3faSDmitry Karpeev x[0][1] = PetscMin(xstart + maxwidth + overlap,M); 1722eec7e3faSDmitry Karpeev /* Horizontal domain limits without an overlap. */ 1723eec7e3faSDmitry Karpeev x[1][0] = xstart; 1724eec7e3faSDmitry Karpeev x[1][1] = PetscMin(xstart+maxwidth,M); 17252bbb417fSDmitry Karpeev /* 17262bbb417fSDmitry Karpeev Determine whether this domain intersects this processor's ownership range of pc->pmat. 17272bbb417fSDmitry Karpeev Do this twice: first for the domains with overlaps, and once without. 17282bbb417fSDmitry Karpeev During the first pass create the subcommunicators, and use them on the second pass as well. 17292bbb417fSDmitry Karpeev */ 17302bbb417fSDmitry Karpeev for (q = 0; q < 2; ++q) { 1731cc96b2e5SBarry Smith PetscBool split = PETSC_FALSE; 17322bbb417fSDmitry Karpeev /* 17332bbb417fSDmitry Karpeev domain limits, (xleft, xright) and (ylow, yheigh) are adjusted 17342bbb417fSDmitry Karpeev according to whether the domain with an overlap or without is considered. 17352bbb417fSDmitry Karpeev */ 17362bbb417fSDmitry Karpeev xleft = x[q][0]; xright = x[q][1]; 17372bbb417fSDmitry Karpeev ylow = y[q][0]; yhigh = y[q][1]; 1738c3518bceSDmitry Karpeev PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,(&xleft_loc),(&ylow_loc),(&xright_loc),(&yhigh_loc),(&nidx)); 1739af538404SDmitry Karpeev nidx *= dof; 1740eec7e3faSDmitry Karpeev n[q] = nidx; 17412bbb417fSDmitry Karpeev /* 1742eec7e3faSDmitry Karpeev Based on the counted number of indices in the local domain *with an overlap*, 1743af538404SDmitry Karpeev construct a subcommunicator of all the processors supporting this domain. 1744eec7e3faSDmitry Karpeev Observe that a domain with an overlap might have nontrivial local support, 1745eec7e3faSDmitry Karpeev while the domain without an overlap might not. Hence, the decision to participate 1746eec7e3faSDmitry Karpeev in the subcommunicator must be based on the domain with an overlap. 17472bbb417fSDmitry Karpeev */ 17482bbb417fSDmitry Karpeev if (q == 0) { 17492fa5cd67SKarl Rupp if (nidx) color = 1; 17502fa5cd67SKarl Rupp else color = MPI_UNDEFINED; 17512bbb417fSDmitry Karpeev ierr = MPI_Comm_split(comm, color, rank, &subcomm);CHKERRQ(ierr); 1752cc96b2e5SBarry Smith split = PETSC_TRUE; 17532bbb417fSDmitry Karpeev } 1754af538404SDmitry Karpeev /* 1755eec7e3faSDmitry Karpeev Proceed only if the number of local indices *with an overlap* is nonzero. 1756af538404SDmitry Karpeev */ 1757eec7e3faSDmitry Karpeev if (n[0]) { 17582fa5cd67SKarl Rupp if (q == 0) xis = is; 1759af538404SDmitry Karpeev if (q == 1) { 1760af538404SDmitry Karpeev /* 1761eec7e3faSDmitry Karpeev The IS for the no-overlap subdomain shares a communicator with the overlapping domain. 1762af538404SDmitry Karpeev Moreover, if the overlap is zero, the two ISs are identical. 1763af538404SDmitry Karpeev */ 1764b1a0a8a3SJed Brown if (overlap == 0) { 1765eec7e3faSDmitry Karpeev (*is_local)[s] = (*is)[s]; 17662bbb417fSDmitry Karpeev ierr = PetscObjectReference((PetscObject)(*is)[s]);CHKERRQ(ierr); 17672bbb417fSDmitry Karpeev continue; 1768d34fcf5fSBarry Smith } else { 17696a4f0f73SDmitry Karpeev xis = is_local; 1770eec7e3faSDmitry Karpeev subcomm = ((PetscObject)(*is)[s])->comm; 17712bbb417fSDmitry Karpeev } 1772af538404SDmitry Karpeev } /* if (q == 1) */ 17730298fd71SBarry Smith idx = NULL; 1774785e854fSJed Brown ierr = PetscMalloc1(nidx,&idx);CHKERRQ(ierr); 1775eec7e3faSDmitry Karpeev if (nidx) { 17762bbb417fSDmitry Karpeev k = 0; 17772bbb417fSDmitry Karpeev for (jj=ylow_loc; jj<yhigh_loc; ++jj) { 1778af538404SDmitry Karpeev PetscInt x0 = (jj==ylow_loc) ? xleft_loc : xleft; 1779af538404SDmitry Karpeev PetscInt x1 = (jj==yhigh_loc-1) ? xright_loc : xright; 1780af538404SDmitry Karpeev kk = dof*(M*jj + x0); 17812bbb417fSDmitry Karpeev for (ii=x0; ii<x1; ++ii) { 17822bbb417fSDmitry Karpeev for (d = 0; d < dof; ++d) { 17832bbb417fSDmitry Karpeev idx[k++] = kk++; 1784b1a0a8a3SJed Brown } 1785b1a0a8a3SJed Brown } 1786b1a0a8a3SJed Brown } 1787eec7e3faSDmitry Karpeev } 17886a4f0f73SDmitry Karpeev ierr = ISCreateGeneral(subcomm,nidx,idx,PETSC_OWN_POINTER,(*xis)+s);CHKERRQ(ierr); 1789cc96b2e5SBarry Smith if (split) { 1790cc96b2e5SBarry Smith ierr = MPI_Comm_free(&subcomm);CHKERRQ(ierr); 1791cc96b2e5SBarry Smith } 1792eec7e3faSDmitry Karpeev }/* if (n[0]) */ 17932bbb417fSDmitry Karpeev }/* for (q = 0; q < 2; ++q) */ 17942fa5cd67SKarl Rupp if (n[0]) ++s; 1795af538404SDmitry Karpeev xstart += maxwidth; 1796af538404SDmitry Karpeev } /* for (i = 0; i < Mdomains; ++i) */ 17972bbb417fSDmitry Karpeev ystart += maxheight; 1798af538404SDmitry Karpeev } /* for (j = 0; j < Ndomains; ++j) */ 1799b1a0a8a3SJed Brown PetscFunctionReturn(0); 1800b1a0a8a3SJed Brown } 1801b1a0a8a3SJed Brown 1802b1a0a8a3SJed Brown /*@C 180306b43e7eSDmitry Karpeev PCGASMGetSubdomains - Gets the subdomains supported on this processor 180406b43e7eSDmitry Karpeev for the additive Schwarz preconditioner. 1805b1a0a8a3SJed Brown 1806b1a0a8a3SJed Brown Not Collective 1807b1a0a8a3SJed Brown 1808b1a0a8a3SJed Brown Input Parameter: 1809b1a0a8a3SJed Brown . pc - the preconditioner context 1810b1a0a8a3SJed Brown 1811b1a0a8a3SJed Brown Output Parameters: 1812b1a0a8a3SJed Brown + n - the number of subdomains for this processor (default value = 1) 18130298fd71SBarry Smith . iis - the index sets that define the inner subdomains (without overlap) supported on this processor (can be NULL) 18140298fd71SBarry Smith - ois - the index sets that define the outer subdomains (with overlap) supported on this processor (can be NULL) 1815b1a0a8a3SJed Brown 1816b1a0a8a3SJed Brown 1817b1a0a8a3SJed Brown Notes: 18186a4f0f73SDmitry Karpeev The user is responsible for destroying the ISs and freeing the returned arrays. 1819b1a0a8a3SJed Brown The IS numbering is in the parallel, global numbering of the vector. 1820b1a0a8a3SJed Brown 1821b1a0a8a3SJed Brown Level: advanced 1822b1a0a8a3SJed Brown 182306b43e7eSDmitry Karpeev .keywords: PC, GASM, get, subdomains, additive Schwarz 1824b1a0a8a3SJed Brown 18258f3b4b4dSDmitry Karpeev .seealso: PCGASMSetOverlap(), PCGASMGetSubKSP(), PCGASMCreateSubdomains2D(), 18268f3b4b4dSDmitry Karpeev PCGASMSetSubdomains(), PCGASMGetSubmatrices() 1827b1a0a8a3SJed Brown @*/ 18286a4f0f73SDmitry Karpeev PetscErrorCode PCGASMGetSubdomains(PC pc,PetscInt *n,IS *iis[],IS *ois[]) 1829b1a0a8a3SJed Brown { 1830f746d493SDmitry Karpeev PC_GASM *osm; 1831b1a0a8a3SJed Brown PetscErrorCode ierr; 1832b1a0a8a3SJed Brown PetscBool match; 18336a4f0f73SDmitry Karpeev PetscInt i; 18345fd66863SKarl Rupp 1835b1a0a8a3SJed Brown PetscFunctionBegin; 1836b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1837251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match);CHKERRQ(ierr); 1838ce94432eSBarry Smith if (!match) SETERRQ2(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Incorrect object type: expected %s, got %s instead", PCGASM, ((PetscObject)pc)->type_name); 1839f746d493SDmitry Karpeev osm = (PC_GASM*)pc->data; 1840ab3e923fSDmitry Karpeev if (n) *n = osm->n; 18416a4f0f73SDmitry Karpeev if (iis) { 1842785e854fSJed Brown ierr = PetscMalloc1(osm->n, iis);CHKERRQ(ierr); 18436a4f0f73SDmitry Karpeev } 18446a4f0f73SDmitry Karpeev if (ois) { 1845785e854fSJed Brown ierr = PetscMalloc1(osm->n, ois);CHKERRQ(ierr); 18466a4f0f73SDmitry Karpeev } 18476a4f0f73SDmitry Karpeev if (iis || ois) { 18486a4f0f73SDmitry Karpeev for (i = 0; i < osm->n; ++i) { 18496a4f0f73SDmitry Karpeev if (iis) (*iis)[i] = osm->iis[i]; 18506a4f0f73SDmitry Karpeev if (ois) (*ois)[i] = osm->ois[i]; 18516a4f0f73SDmitry Karpeev } 1852b1a0a8a3SJed Brown } 1853b1a0a8a3SJed Brown PetscFunctionReturn(0); 1854b1a0a8a3SJed Brown } 1855b1a0a8a3SJed Brown 1856b1a0a8a3SJed Brown /*@C 185706b43e7eSDmitry Karpeev PCGASMGetSubmatrices - Gets the local submatrices (for this processor 1858b1a0a8a3SJed Brown only) for the additive Schwarz preconditioner. 1859b1a0a8a3SJed Brown 1860b1a0a8a3SJed Brown Not Collective 1861b1a0a8a3SJed Brown 1862b1a0a8a3SJed Brown Input Parameter: 1863b1a0a8a3SJed Brown . pc - the preconditioner context 1864b1a0a8a3SJed Brown 1865b1a0a8a3SJed Brown Output Parameters: 1866b1a0a8a3SJed Brown + n - the number of matrices for this processor (default value = 1) 1867b1a0a8a3SJed Brown - mat - the matrices 1868b1a0a8a3SJed Brown 186995452b02SPatrick Sanan Notes: 187095452b02SPatrick Sanan matrices returned by this routine have the same communicators as the index sets (IS) 18718f3b4b4dSDmitry Karpeev used to define subdomains in PCGASMSetSubdomains() 1872b1a0a8a3SJed Brown Level: advanced 1873b1a0a8a3SJed Brown 1874f746d493SDmitry Karpeev .keywords: PC, GASM, set, local, subdomains, additive Schwarz, block Jacobi 1875b1a0a8a3SJed Brown 18768f3b4b4dSDmitry Karpeev .seealso: PCGASMSetOverlap(), PCGASMGetSubKSP(), 187706b43e7eSDmitry Karpeev PCGASMCreateSubdomains2D(), PCGASMSetSubdomains(), PCGASMGetSubdomains() 1878b1a0a8a3SJed Brown @*/ 187906b43e7eSDmitry Karpeev PetscErrorCode PCGASMGetSubmatrices(PC pc,PetscInt *n,Mat *mat[]) 1880b1a0a8a3SJed Brown { 1881f746d493SDmitry Karpeev PC_GASM *osm; 1882b1a0a8a3SJed Brown PetscErrorCode ierr; 1883b1a0a8a3SJed Brown PetscBool match; 1884b1a0a8a3SJed Brown 1885b1a0a8a3SJed Brown PetscFunctionBegin; 1886b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1887b1a0a8a3SJed Brown PetscValidIntPointer(n,2); 1888b1a0a8a3SJed Brown if (mat) PetscValidPointer(mat,3); 1889ce94432eSBarry Smith if (!pc->setupcalled) SETERRQ(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call after KSPSetUP() or PCSetUp()."); 1890251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match);CHKERRQ(ierr); 1891ce94432eSBarry Smith if (!match) SETERRQ2(PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Expected %s, got %s instead", PCGASM, ((PetscObject)pc)->type_name); 1892f746d493SDmitry Karpeev osm = (PC_GASM*)pc->data; 1893ab3e923fSDmitry Karpeev if (n) *n = osm->n; 1894b1a0a8a3SJed Brown if (mat) *mat = osm->pmat; 1895b1a0a8a3SJed Brown PetscFunctionReturn(0); 1896b1a0a8a3SJed Brown } 1897d709ea83SDmitry Karpeev 1898d709ea83SDmitry Karpeev /*@ 18998f3b4b4dSDmitry Karpeev PCGASMSetUseDMSubdomains - Indicates whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible. 1900d709ea83SDmitry Karpeev Logically Collective 1901d709ea83SDmitry Karpeev 1902d709ea83SDmitry Karpeev Input Parameter: 1903d709ea83SDmitry Karpeev + pc - the preconditioner 1904d709ea83SDmitry Karpeev - flg - boolean indicating whether to use subdomains defined by the DM 1905d709ea83SDmitry Karpeev 1906d709ea83SDmitry Karpeev Options Database Key: 19078f3b4b4dSDmitry Karpeev . -pc_gasm_dm_subdomains -pc_gasm_overlap -pc_gasm_total_subdomains 1908d709ea83SDmitry Karpeev 1909d709ea83SDmitry Karpeev Level: intermediate 1910d709ea83SDmitry Karpeev 1911d709ea83SDmitry Karpeev Notes: 19128f3b4b4dSDmitry Karpeev PCGASMSetSubdomains(), PCGASMSetTotalSubdomains() or PCGASMSetOverlap() take precedence over PCGASMSetUseDMSubdomains(), 19138f3b4b4dSDmitry Karpeev so setting PCGASMSetSubdomains() with nontrivial subdomain ISs or any of PCGASMSetTotalSubdomains() and PCGASMSetOverlap() 19148f3b4b4dSDmitry Karpeev automatically turns the latter off. 1915d709ea83SDmitry Karpeev 1916d709ea83SDmitry Karpeev .keywords: PC, ASM, DM, set, subdomains, additive Schwarz 1917d709ea83SDmitry Karpeev 19188f3b4b4dSDmitry Karpeev .seealso: PCGASMGetUseDMSubdomains(), PCGASMSetSubdomains(), PCGASMSetOverlap() 1919d709ea83SDmitry Karpeev PCGASMCreateSubdomains2D() 1920d709ea83SDmitry Karpeev @*/ 19218f3b4b4dSDmitry Karpeev PetscErrorCode PCGASMSetUseDMSubdomains(PC pc,PetscBool flg) 1922d709ea83SDmitry Karpeev { 1923d709ea83SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 1924d709ea83SDmitry Karpeev PetscErrorCode ierr; 1925d709ea83SDmitry Karpeev PetscBool match; 1926d709ea83SDmitry Karpeev 1927d709ea83SDmitry Karpeev PetscFunctionBegin; 1928d709ea83SDmitry Karpeev PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1929d709ea83SDmitry Karpeev PetscValidLogicalCollectiveBool(pc,flg,2); 1930d709ea83SDmitry Karpeev if (pc->setupcalled) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Not for a setup PC."); 1931d709ea83SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match);CHKERRQ(ierr); 1932d709ea83SDmitry Karpeev if (match) { 19338f3b4b4dSDmitry Karpeev if (!osm->user_subdomains && osm->N == PETSC_DETERMINE && osm->overlap < 0) { 1934d709ea83SDmitry Karpeev osm->dm_subdomains = flg; 1935d709ea83SDmitry Karpeev } 19368f3b4b4dSDmitry Karpeev } 1937d709ea83SDmitry Karpeev PetscFunctionReturn(0); 1938d709ea83SDmitry Karpeev } 1939d709ea83SDmitry Karpeev 1940d709ea83SDmitry Karpeev /*@ 19418f3b4b4dSDmitry Karpeev PCGASMGetUseDMSubdomains - Returns flag indicating whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible. 1942d709ea83SDmitry Karpeev Not Collective 1943d709ea83SDmitry Karpeev 1944d709ea83SDmitry Karpeev Input Parameter: 1945d709ea83SDmitry Karpeev . pc - the preconditioner 1946d709ea83SDmitry Karpeev 1947d709ea83SDmitry Karpeev Output Parameter: 1948d709ea83SDmitry Karpeev . flg - boolean indicating whether to use subdomains defined by the DM 1949d709ea83SDmitry Karpeev 1950d709ea83SDmitry Karpeev Level: intermediate 1951d709ea83SDmitry Karpeev 1952d709ea83SDmitry Karpeev .keywords: PC, ASM, DM, set, subdomains, additive Schwarz 1953d709ea83SDmitry Karpeev 19548f3b4b4dSDmitry Karpeev .seealso: PCGASMSetUseDMSubdomains(), PCGASMSetOverlap() 1955d709ea83SDmitry Karpeev PCGASMCreateSubdomains2D() 1956d709ea83SDmitry Karpeev @*/ 19578f3b4b4dSDmitry Karpeev PetscErrorCode PCGASMGetUseDMSubdomains(PC pc,PetscBool* flg) 1958d709ea83SDmitry Karpeev { 1959d709ea83SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 1960d709ea83SDmitry Karpeev PetscErrorCode ierr; 1961d709ea83SDmitry Karpeev PetscBool match; 1962d709ea83SDmitry Karpeev 1963d709ea83SDmitry Karpeev PetscFunctionBegin; 1964d709ea83SDmitry Karpeev PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1965d709ea83SDmitry Karpeev PetscValidPointer(flg,2); 1966d709ea83SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match);CHKERRQ(ierr); 1967d709ea83SDmitry Karpeev if (match) { 1968d709ea83SDmitry Karpeev if (flg) *flg = osm->dm_subdomains; 1969d709ea83SDmitry Karpeev } 1970d709ea83SDmitry Karpeev PetscFunctionReturn(0); 1971d709ea83SDmitry Karpeev } 1972