1b1a0a8a3SJed Brown /* 2f746d493SDmitry Karpeev This file defines an "generalized" additive Schwarz preconditioner for any Mat implementation. 36a4f0f73SDmitry Karpeev In this version each processor may intersect multiple subdomains and any subdomain may 46a4f0f73SDmitry Karpeev intersect multiple processors. Intersections of subdomains with processors are called *local 56a4f0f73SDmitry Karpeev subdomains*. 6b1a0a8a3SJed Brown 76a4f0f73SDmitry Karpeev N - total number of local subdomains on all processors (set in PCGASMSetTotalSubdomains() or calculated in PCSetUp_GASM()) 86a4f0f73SDmitry Karpeev n - actual number of local subdomains on this processor (set in PCGASMSetSubdomains() or calculated in PCGASMSetTotalSubdomains()) 96a4f0f73SDmitry Karpeev nmax - maximum number of local subdomains per processor (calculated in PCGASMSetTotalSubdomains() or in PCSetUp_GASM()) 10b1a0a8a3SJed Brown */ 11b45d2f2cSJed Brown #include <petsc-private/pcimpl.h> /*I "petscpc.h" I*/ 12b1a0a8a3SJed Brown 13b1a0a8a3SJed Brown typedef struct { 14f746d493SDmitry Karpeev PetscInt N,n,nmax; 15b1a0a8a3SJed Brown PetscInt overlap; /* overlap requested by user */ 16b1a0a8a3SJed Brown KSP *ksp; /* linear solvers for each block */ 17f746d493SDmitry Karpeev Vec gx,gy; /* Merged work vectors */ 18f746d493SDmitry Karpeev Vec *x,*y; /* Split work vectors; storage aliases pieces of storage of the above merged vectors. */ 196a4f0f73SDmitry Karpeev VecScatter gorestriction; /* merged restriction to disjoint union of outer subdomains */ 206a4f0f73SDmitry Karpeev VecScatter girestriction; /* merged restriction to disjoint union of inner subdomains */ 216a4f0f73SDmitry Karpeev IS *ois; /* index sets that define the outer (conceptually, overlapping) subdomains */ 226a4f0f73SDmitry Karpeev IS *iis; /* index sets that define the inner (conceptually, nonoverlapping) subdomains */ 23f746d493SDmitry Karpeev Mat *pmat; /* subdomain block matrices */ 24f746d493SDmitry Karpeev PCGASMType type; /* use reduced interpolation, restriction or both */ 256a4f0f73SDmitry Karpeev PetscBool create_local; /* whether the autocreated subdomains are local or not. */ 26b1a0a8a3SJed Brown PetscBool type_set; /* if user set this value (so won't change it for symmetric problems) */ 276a4f0f73SDmitry Karpeev PetscBool same_subdomain_solvers; /* flag indicating whether all local solvers are same */ 28b1a0a8a3SJed Brown PetscBool sort_indices; /* flag to sort subdomain indices */ 29f746d493SDmitry Karpeev } PC_GASM; 30b1a0a8a3SJed Brown 31b1a0a8a3SJed Brown #undef __FUNCT__ 3206b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMSubdomainView_Private" 3306b43e7eSDmitry Karpeev static PetscErrorCode PCGASMSubdomainView_Private(PC pc, PetscInt i, PetscViewer viewer) 34af538404SDmitry Karpeev { 35af538404SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 3606b43e7eSDmitry Karpeev PetscInt j,nidx; 37af538404SDmitry Karpeev const PetscInt *idx; 3806b43e7eSDmitry Karpeev PetscViewer sviewer; 39af538404SDmitry Karpeev char *cidx; 40af538404SDmitry Karpeev PetscErrorCode ierr; 41af538404SDmitry Karpeev 42af538404SDmitry Karpeev PetscFunctionBegin; 4306b43e7eSDmitry Karpeev if (i < 0 || i > osm->n) SETERRQ2(((PetscObject)viewer)->comm, PETSC_ERR_ARG_WRONG, "Invalid subdomain %D: must nonnegative and less than %D", i, osm->n); 444bde246dSDmitry Karpeev /* Inner subdomains. */ 454bde246dSDmitry Karpeev ierr = ISGetLocalSize(osm->iis[i], &nidx);CHKERRQ(ierr); 464bde246dSDmitry Karpeev /* 474bde246dSDmitry Karpeev No more than 15 characters per index plus a space. 484bde246dSDmitry Karpeev PetscViewerStringSPrintf requires a string of size at least 2, so use (nidx+1) instead of nidx, 494bde246dSDmitry Karpeev in case nidx == 0. That will take care of the space for the trailing '\0' as well. 504bde246dSDmitry Karpeev For nidx == 0, the whole string 16 '\0'. 514bde246dSDmitry Karpeev */ 524bde246dSDmitry Karpeev ierr = PetscMalloc(sizeof(char)*(16*(nidx+1)+1), &cidx);CHKERRQ(ierr); 534bde246dSDmitry Karpeev ierr = PetscViewerStringOpen(PETSC_COMM_SELF, cidx, 16*(nidx+1)+1, &sviewer);CHKERRQ(ierr); 544bde246dSDmitry Karpeev ierr = ISGetIndices(osm->iis[i], &idx);CHKERRQ(ierr); 554bde246dSDmitry Karpeev for (j = 0; j < nidx; ++j) { 564bde246dSDmitry Karpeev ierr = PetscViewerStringSPrintf(sviewer, "%D ", idx[j]);CHKERRQ(ierr); 574bde246dSDmitry Karpeev } 584bde246dSDmitry Karpeev ierr = ISRestoreIndices(osm->iis[i],&idx);CHKERRQ(ierr); 594bde246dSDmitry Karpeev ierr = PetscViewerDestroy(&sviewer);CHKERRQ(ierr); 604bde246dSDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer, "Inner subdomain:\n");CHKERRQ(ierr); 614bde246dSDmitry Karpeev ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 624bde246dSDmitry Karpeev ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr); 634bde246dSDmitry Karpeev ierr = PetscViewerASCIISynchronizedPrintf(viewer, "%s", cidx);CHKERRQ(ierr); 644bde246dSDmitry Karpeev ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 654bde246dSDmitry Karpeev ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr); 664bde246dSDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer, "\n");CHKERRQ(ierr); 674bde246dSDmitry Karpeev ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 684bde246dSDmitry Karpeev ierr = PetscFree(cidx);CHKERRQ(ierr); 694bde246dSDmitry Karpeev /* Outer subdomains. */ 706a4f0f73SDmitry Karpeev ierr = ISGetLocalSize(osm->ois[i], &nidx);CHKERRQ(ierr); 71eec7e3faSDmitry Karpeev /* 72eec7e3faSDmitry Karpeev No more than 15 characters per index plus a space. 73eec7e3faSDmitry Karpeev PetscViewerStringSPrintf requires a string of size at least 2, so use (nidx+1) instead of nidx, 74eec7e3faSDmitry Karpeev in case nidx == 0. That will take care of the space for the trailing '\0' as well. 75eec7e3faSDmitry Karpeev For nidx == 0, the whole string 16 '\0'. 76eec7e3faSDmitry Karpeev */ 77eec7e3faSDmitry Karpeev ierr = PetscMalloc(sizeof(char)*(16*(nidx+1)+1), &cidx);CHKERRQ(ierr); 7806b43e7eSDmitry Karpeev ierr = PetscViewerStringOpen(PETSC_COMM_SELF, cidx, 16*(nidx+1)+1, &sviewer);CHKERRQ(ierr); 796a4f0f73SDmitry Karpeev ierr = ISGetIndices(osm->ois[i], &idx);CHKERRQ(ierr); 80af538404SDmitry Karpeev for (j = 0; j < nidx; ++j) { 81af538404SDmitry Karpeev ierr = PetscViewerStringSPrintf(sviewer,"%D ", idx[j]);CHKERRQ(ierr); 82af538404SDmitry Karpeev } 836bf464f9SBarry Smith ierr = PetscViewerDestroy(&sviewer);CHKERRQ(ierr); 846a4f0f73SDmitry Karpeev ierr = ISRestoreIndices(osm->ois[i],&idx);CHKERRQ(ierr); 856a4f0f73SDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer, "Outer subdomain:\n");CHKERRQ(ierr); 8606b43e7eSDmitry Karpeev ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 877b23a99aSBarry Smith ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr); 88af538404SDmitry Karpeev ierr = PetscViewerASCIISynchronizedPrintf(viewer, "%s", cidx);CHKERRQ(ierr); 89af538404SDmitry Karpeev ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 907b23a99aSBarry Smith ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr); 91af538404SDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer, "\n");CHKERRQ(ierr); 9206b43e7eSDmitry Karpeev ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 93af538404SDmitry Karpeev ierr = PetscFree(cidx);CHKERRQ(ierr); 944bde246dSDmitry Karpeev 9506b43e7eSDmitry Karpeev PetscFunctionReturn(0); 9606b43e7eSDmitry Karpeev } 9706b43e7eSDmitry Karpeev 9806b43e7eSDmitry Karpeev #undef __FUNCT__ 9906b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMPrintSubdomains" 10006b43e7eSDmitry Karpeev static PetscErrorCode PCGASMPrintSubdomains(PC pc) 10106b43e7eSDmitry Karpeev { 10206b43e7eSDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 10306b43e7eSDmitry Karpeev const char *prefix; 10406b43e7eSDmitry Karpeev char fname[PETSC_MAX_PATH_LEN+1]; 1054bde246dSDmitry Karpeev PetscInt i, l, d, count, gcount, *permutation, *numbering; 10606b43e7eSDmitry Karpeev PetscBool found; 1074bde246dSDmitry Karpeev PetscViewer viewer, sviewer = PETSC_NULL; 10806b43e7eSDmitry Karpeev PetscErrorCode ierr; 10906b43e7eSDmitry Karpeev 11006b43e7eSDmitry Karpeev PetscFunctionBegin; 1114bde246dSDmitry Karpeev ierr = PetscMalloc2(osm->n, PetscInt, &permutation, osm->n, PetscInt, &numbering);CHKERRQ(ierr); 1124bde246dSDmitry Karpeev for (i = 0; i < osm->n; ++i) permutation[i] = i; 1134bde246dSDmitry Karpeev ierr = PetscObjectsGetGlobalNumbering(((PetscObject)pc)->comm, osm->n, (PetscObject*)osm->ois, &gcount, numbering);CHKERRQ(ierr); 1144bde246dSDmitry Karpeev ierr = PetscSortIntWithPermutation(osm->n, numbering, permutation);CHKERRQ(ierr); 11506b43e7eSDmitry Karpeev ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 11606b43e7eSDmitry Karpeev ierr = PetscOptionsGetString(prefix,"-pc_gasm_print_subdomains",fname,PETSC_MAX_PATH_LEN,&found);CHKERRQ(ierr); 11706b43e7eSDmitry Karpeev if (!found) { ierr = PetscStrcpy(fname,"stdout");CHKERRQ(ierr); }; 1184bde246dSDmitry Karpeev ierr = PetscViewerASCIIOpen(((PetscObject)pc)->comm,fname,&viewer);CHKERRQ(ierr); 1194bde246dSDmitry Karpeev /* 1204bde246dSDmitry Karpeev Make sure the viewer has a name. Otherwise this may cause a deadlock or other weird errors when creating a subcomm viewer: 1214bde246dSDmitry Karpeev the subcomm viewer will attempt to inherit the viewer's name, which, if not set, will be constructed collectively on the comm. 1224bde246dSDmitry Karpeev */ 1234bde246dSDmitry Karpeev ierr = PetscObjectName((PetscObject)viewer);CHKERRQ(ierr); 1244bde246dSDmitry Karpeev l = 0; 1254bde246dSDmitry Karpeev for (count = 0; count < gcount; ++count) { 1264bde246dSDmitry Karpeev /* Now let subdomains go one at a time in the global numbering order and print their subdomain/solver info. */ 1274bde246dSDmitry Karpeev if (l<osm->n) { 1284bde246dSDmitry Karpeev d = permutation[l]; /* d is the local number of the l-th smallest (in the global ordering) among the locally supported subdomains */ 1294bde246dSDmitry Karpeev if (numbering[d] == count) { 1304bde246dSDmitry Karpeev ierr = PetscViewerGetSubcomm(viewer,((PetscObject)osm->ois[d])->comm, &sviewer);CHKERRQ(ierr); 1314bde246dSDmitry Karpeev ierr = PCGASMSubdomainView_Private(pc,d,sviewer);CHKERRQ(ierr); 1324bde246dSDmitry Karpeev ierr = PetscViewerRestoreSubcomm(viewer,((PetscObject)osm->ois[d])->comm, &sviewer);CHKERRQ(ierr); 1334bde246dSDmitry Karpeev ++l; 134af538404SDmitry Karpeev } 1354bde246dSDmitry Karpeev } 1364bde246dSDmitry Karpeev ierr = MPI_Barrier(((PetscObject)pc)->comm);CHKERRQ(ierr); 1374bde246dSDmitry Karpeev } 1384bde246dSDmitry Karpeev ierr = PetscViewerDestroy(&viewer);CHKERRQ(ierr); 1397105d80bSDmitry Karpeev ierr = PetscFree2(permutation,numbering);CHKERRQ(ierr); 140af538404SDmitry Karpeev PetscFunctionReturn(0); 141af538404SDmitry Karpeev } 142af538404SDmitry Karpeev 143af538404SDmitry Karpeev 144af538404SDmitry Karpeev #undef __FUNCT__ 145f746d493SDmitry Karpeev #define __FUNCT__ "PCView_GASM" 146f746d493SDmitry Karpeev static PetscErrorCode PCView_GASM(PC pc,PetscViewer viewer) 147b1a0a8a3SJed Brown { 148f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 149af538404SDmitry Karpeev const char *prefix; 150b1a0a8a3SJed Brown PetscErrorCode ierr; 151af538404SDmitry Karpeev PetscMPIInt rank, size; 152b1a0a8a3SJed Brown PetscInt i,bsz; 15306b43e7eSDmitry Karpeev PetscBool iascii,view_subdomains=PETSC_FALSE; 154b1a0a8a3SJed Brown PetscViewer sviewer; 15506b43e7eSDmitry Karpeev PetscInt count, l, gcount, *numbering, *permutation; 1566a4f0f73SDmitry Karpeev char overlap[256] = "user-defined overlap"; 1576a4f0f73SDmitry Karpeev char gsubdomains[256] = "unknown total number of subdomains"; 15806b43e7eSDmitry Karpeev char lsubdomains[256] = "unknown number of local subdomains"; 15906b43e7eSDmitry Karpeev char msubdomains[256] = "unknown max number of local subdomains"; 160b1a0a8a3SJed Brown PetscFunctionBegin; 161af538404SDmitry Karpeev ierr = MPI_Comm_size(((PetscObject)pc)->comm, &size);CHKERRQ(ierr); 162af538404SDmitry Karpeev ierr = MPI_Comm_rank(((PetscObject)pc)->comm, &rank);CHKERRQ(ierr); 16306b43e7eSDmitry Karpeev 16406b43e7eSDmitry Karpeev 16506b43e7eSDmitry Karpeev ierr = PetscMalloc2(osm->n, PetscInt, &permutation, osm->n, PetscInt, &numbering);CHKERRQ(ierr); 16606b43e7eSDmitry Karpeev for (i = 0; i < osm->n; ++i) permutation[i] = i; 1676a4f0f73SDmitry Karpeev ierr = PetscObjectsGetGlobalNumbering(((PetscObject)pc)->comm, osm->n, (PetscObject*)osm->ois, &gcount, numbering);CHKERRQ(ierr); 16806b43e7eSDmitry Karpeev ierr = PetscSortIntWithPermutation(osm->n, numbering, permutation);CHKERRQ(ierr); 16906b43e7eSDmitry Karpeev 17006b43e7eSDmitry Karpeev if (osm->overlap >= 0) { 17106b43e7eSDmitry Karpeev ierr = PetscSNPrintf(overlap,sizeof(overlap),"requested amount of overlap = %D",osm->overlap);CHKERRQ(ierr); 17206b43e7eSDmitry Karpeev } 1736a4f0f73SDmitry Karpeev ierr = PetscSNPrintf(gsubdomains, sizeof(gsubdomains), "total number of subdomains = %D",gcount);CHKERRQ(ierr); 17406b43e7eSDmitry Karpeev if (osm->N > 0) { 17506b43e7eSDmitry Karpeev ierr = PetscSNPrintf(lsubdomains, sizeof(gsubdomains), "number of local subdomains = %D",osm->N);CHKERRQ(ierr); 17606b43e7eSDmitry Karpeev } 17706b43e7eSDmitry Karpeev if (osm->nmax > 0) { 17806b43e7eSDmitry Karpeev ierr = PetscSNPrintf(msubdomains,sizeof(msubdomains),"max number of local subdomains = %D",osm->nmax);CHKERRQ(ierr); 17906b43e7eSDmitry Karpeev } 18006b43e7eSDmitry Karpeev 181af538404SDmitry Karpeev ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 18206b43e7eSDmitry Karpeev ierr = PetscOptionsGetBool(prefix,"-pc_gasm_view_subdomains",&view_subdomains,PETSC_NULL);CHKERRQ(ierr); 18306b43e7eSDmitry Karpeev 18406b43e7eSDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 185b1a0a8a3SJed Brown if (iascii) { 18606b43e7eSDmitry Karpeev /* 18706b43e7eSDmitry Karpeev Make sure the viewer has a name. Otherwise this may cause a deadlock when creating a subcomm viewer: 18806b43e7eSDmitry Karpeev the subcomm viewer will attempt to inherit the viewer's name, which, if not set, will be constructed 18906b43e7eSDmitry Karpeev collectively on the comm. 19006b43e7eSDmitry Karpeev */ 19106b43e7eSDmitry Karpeev ierr = PetscObjectName((PetscObject)viewer);CHKERRQ(ierr); 19206b43e7eSDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer,"Generalized additive Schwarz:\n");CHKERRQ(ierr); 19306b43e7eSDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer,"Restriction/interpolation type: %s\n",PCGASMTypes[osm->type]);CHKERRQ(ierr); 19406b43e7eSDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer,"%s\n",overlap);CHKERRQ(ierr); 19506b43e7eSDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer,"%s\n",gsubdomains);CHKERRQ(ierr); 19606b43e7eSDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer,"%s\n",lsubdomains);CHKERRQ(ierr); 19706b43e7eSDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer,"%s\n",msubdomains);CHKERRQ(ierr); 1987b23a99aSBarry Smith ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr); 199eec7e3faSDmitry Karpeev ierr = PetscViewerASCIISynchronizedPrintf(viewer,"[%d:%d] number of locally-supported subdomains = %D\n",(int)rank,(int)size,osm->n);CHKERRQ(ierr); 200af538404SDmitry Karpeev ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 2017b23a99aSBarry Smith ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr); 2026a4f0f73SDmitry Karpeev /* Cannot take advantage of osm->same_subdomain_solvers without a global numbering of subdomains. */ 20306b43e7eSDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer,"Subdomain solver info is as follows:\n");CHKERRQ(ierr); 204b1a0a8a3SJed Brown ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 205b1a0a8a3SJed Brown ierr = PetscViewerASCIIPrintf(viewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr); 20606b43e7eSDmitry Karpeev /* Make sure that everybody waits for the banner to be printed. */ 20706b43e7eSDmitry Karpeev ierr = MPI_Barrier(((PetscObject)viewer)->comm);CHKERRQ(ierr); 2084bde246dSDmitry Karpeev /* Now let subdomains go one at a time in the global numbering order and print their subdomain/solver info. */ 20906b43e7eSDmitry Karpeev l = 0; 21006b43e7eSDmitry Karpeev for (count = 0; count < gcount; ++count) { 21106b43e7eSDmitry Karpeev PetscMPIInt srank, ssize; 21206b43e7eSDmitry Karpeev if (l<osm->n) { 21306b43e7eSDmitry Karpeev PetscInt d = permutation[l]; /* d is the local number of the l-th smallest (in the global ordering) among the locally supported subdomains */ 21406b43e7eSDmitry Karpeev if (numbering[d] == count) { 2156a4f0f73SDmitry Karpeev ierr = MPI_Comm_size(((PetscObject)osm->ois[d])->comm, &ssize);CHKERRQ(ierr); 2166a4f0f73SDmitry Karpeev ierr = MPI_Comm_rank(((PetscObject)osm->ois[d])->comm, &srank);CHKERRQ(ierr); 2176a4f0f73SDmitry Karpeev ierr = PetscViewerGetSubcomm(viewer,((PetscObject)osm->ois[d])->comm, &sviewer);CHKERRQ(ierr); 2186a4f0f73SDmitry Karpeev ierr = ISGetLocalSize(osm->ois[d],&bsz);CHKERRQ(ierr); 2197b23a99aSBarry Smith ierr = PetscViewerASCIISynchronizedAllow(sviewer,PETSC_TRUE);CHKERRQ(ierr); 22006b43e7eSDmitry Karpeev ierr = PetscViewerASCIISynchronizedPrintf(sviewer,"[%D:%D] (subcomm [%D:%D]) local subdomain number %D, local size = %D\n",(int)rank,(int)size,(int)srank,(int)ssize,d,bsz);CHKERRQ(ierr); 2216a4f0f73SDmitry Karpeev ierr = PetscViewerFlush(sviewer);CHKERRQ(ierr); 2226a4f0f73SDmitry Karpeev ierr = PetscViewerASCIISynchronizedAllow(sviewer,PETSC_FALSE);CHKERRQ(ierr); 22306b43e7eSDmitry Karpeev if (view_subdomains) { 22406b43e7eSDmitry Karpeev ierr = PCGASMSubdomainView_Private(pc,d,sviewer);CHKERRQ(ierr); 22506b43e7eSDmitry Karpeev } 2266a4f0f73SDmitry Karpeev if (!pc->setupcalled) { 2276a4f0f73SDmitry Karpeev PetscViewerASCIIPrintf(sviewer, "Solver not set up yet: PCSetUp() not yet called\n");CHKERRQ(ierr); 2286a4f0f73SDmitry Karpeev } 2296a4f0f73SDmitry Karpeev else { 23006b43e7eSDmitry Karpeev ierr = KSPView(osm->ksp[d],sviewer);CHKERRQ(ierr); 2316a4f0f73SDmitry Karpeev } 23206b43e7eSDmitry Karpeev ierr = PetscViewerASCIIPrintf(sviewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr); 2337b23a99aSBarry Smith ierr = PetscViewerFlush(sviewer);CHKERRQ(ierr); 2346a4f0f73SDmitry Karpeev ierr = PetscViewerRestoreSubcomm(viewer,((PetscObject)osm->ois[d])->comm, &sviewer);CHKERRQ(ierr); 23506b43e7eSDmitry Karpeev ++l; 236b1a0a8a3SJed Brown } 23706b43e7eSDmitry Karpeev } 23806b43e7eSDmitry Karpeev ierr = MPI_Barrier(((PetscObject)pc)->comm);CHKERRQ(ierr); 239b1a0a8a3SJed Brown } 240b1a0a8a3SJed Brown ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 241b1a0a8a3SJed Brown } else { 242f746d493SDmitry Karpeev SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCGASM",((PetscObject)viewer)->type_name); 243b1a0a8a3SJed Brown } 24406b43e7eSDmitry Karpeev ierr = PetscFree2(permutation,numbering);CHKERRQ(ierr); 245b1a0a8a3SJed Brown PetscFunctionReturn(0); 246b1a0a8a3SJed Brown } 247b1a0a8a3SJed Brown 248af538404SDmitry Karpeev 249af538404SDmitry Karpeev 250af538404SDmitry Karpeev 251af538404SDmitry Karpeev 252b1a0a8a3SJed Brown #undef __FUNCT__ 253f746d493SDmitry Karpeev #define __FUNCT__ "PCSetUp_GASM" 254f746d493SDmitry Karpeev static PetscErrorCode PCSetUp_GASM(PC pc) 255b1a0a8a3SJed Brown { 256f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 257b1a0a8a3SJed Brown PetscErrorCode ierr; 258b1a0a8a3SJed Brown PetscBool symset,flg; 25988389755SJed Brown PetscInt i; 260c10bc1a1SDmitry Karpeev PetscMPIInt rank, size; 261b1a0a8a3SJed Brown MatReuse scall = MAT_REUSE_MATRIX; 262b1a0a8a3SJed Brown KSP ksp; 263b1a0a8a3SJed Brown PC subpc; 264b1a0a8a3SJed Brown const char *prefix,*pprefix; 265f746d493SDmitry Karpeev Vec x,y; 2666a4f0f73SDmitry Karpeev PetscInt oni; /* Number of indices in the i-th local outer subdomain. */ 2676a4f0f73SDmitry Karpeev const PetscInt *oidxi; /* Indices from the i-th subdomain local outer subdomain. */ 2686a4f0f73SDmitry Karpeev PetscInt on; /* Number of indices in the disjoint union of local outer subdomains. */ 2696a4f0f73SDmitry Karpeev PetscInt *oidx; /* Indices in the disjoint union of local outer subdomains. */ 2706a4f0f73SDmitry Karpeev IS gois; /* Disjoint union the global indices of outer subdomains. */ 2716a4f0f73SDmitry Karpeev IS goid; /* Identity IS of the size of the disjoint union of outer subdomains. */ 272f746d493SDmitry Karpeev PetscScalar *gxarray, *gyarray; 2736a4f0f73SDmitry Karpeev PetscInt gofirst; /* Start of locally-owned indices in the vectors -- osm->gx,osm->gy -- 2746a4f0f73SDmitry Karpeev over the disjoint union of outer subdomains. */ 275a35b7b57SDmitry Karpeev DM *subdomain_dm = PETSC_NULL; 276b1a0a8a3SJed Brown 277b1a0a8a3SJed Brown PetscFunctionBegin; 278c10bc1a1SDmitry Karpeev ierr = MPI_Comm_size(((PetscObject)pc)->comm,&size);CHKERRQ(ierr); 279c10bc1a1SDmitry Karpeev ierr = MPI_Comm_rank(((PetscObject)pc)->comm,&rank);CHKERRQ(ierr); 280b1a0a8a3SJed Brown if (!pc->setupcalled) { 281b1a0a8a3SJed Brown 282b1a0a8a3SJed Brown if (!osm->type_set) { 283b1a0a8a3SJed Brown ierr = MatIsSymmetricKnown(pc->pmat,&symset,&flg);CHKERRQ(ierr); 284f746d493SDmitry Karpeev if (symset && flg) { osm->type = PC_GASM_BASIC; } 285b1a0a8a3SJed Brown } 286b1a0a8a3SJed Brown 28706b43e7eSDmitry Karpeev /* 28806b43e7eSDmitry Karpeev If subdomains have been set, then the local number of subdomains, osm->n, is NOT PETSC_DECIDE and is at least 1. 28906b43e7eSDmitry Karpeev The total number of subdomains, osm->N is not necessarily set, might be PETSC_DECIDE, and then will have to be calculated from osm->n. 29006b43e7eSDmitry Karpeev */ 29106b43e7eSDmitry Karpeev if (osm->n == PETSC_DECIDE) { 29269ca1f37SDmitry Karpeev /* no subdomains given */ 29369ca1f37SDmitry Karpeev /* try pc->dm first */ 294fdc48646SDmitry Karpeev if (pc->dm) { 295fdc48646SDmitry Karpeev char ddm_name[1024]; 296fdc48646SDmitry Karpeev DM ddm; 297fdc48646SDmitry Karpeev PetscBool flg; 298a35b7b57SDmitry Karpeev PetscInt num_subdomains, d; 299a35b7b57SDmitry Karpeev char **subdomain_names; 300a35b7b57SDmitry Karpeev IS *inner_subdomain_is, *outer_subdomain_is; 301fdc48646SDmitry Karpeev /* Allow the user to request a decomposition DM by name */ 302fdc48646SDmitry Karpeev ierr = PetscStrncpy(ddm_name, "", 1024);CHKERRQ(ierr); 303a35b7b57SDmitry Karpeev ierr = PetscOptionsString("-pc_gasm_decomposition","Name of the DM defining the decomposition", "PCSetDM",ddm_name,ddm_name,1024,&flg);CHKERRQ(ierr); 304fdc48646SDmitry Karpeev if (flg) { 30516621825SDmitry Karpeev ierr = DMCreateDomainDecompositionDM(pc->dm, ddm_name, &ddm);CHKERRQ(ierr); 306fdc48646SDmitry Karpeev if (!ddm) { 307fdc48646SDmitry Karpeev SETERRQ1(((PetscObject)pc)->comm, PETSC_ERR_ARG_WRONGSTATE, "Uknown DM decomposition name %s", ddm_name); 308fdc48646SDmitry Karpeev } 309fdc48646SDmitry Karpeev ierr = PetscInfo(pc,"Using decomposition DM defined using options database\n");CHKERRQ(ierr); 310fdc48646SDmitry Karpeev ierr = PCSetDM(pc,ddm);CHKERRQ(ierr); 311fdc48646SDmitry Karpeev } 312a35b7b57SDmitry Karpeev ierr = DMCreateDomainDecomposition(pc->dm, &num_subdomains, &subdomain_names, &inner_subdomain_is, &outer_subdomain_is, &subdomain_dm);CHKERRQ(ierr); 313a35b7b57SDmitry Karpeev if (num_subdomains) { 314a35b7b57SDmitry Karpeev ierr = PCGASMSetSubdomains(pc, num_subdomains, inner_subdomain_is, outer_subdomain_is);CHKERRQ(ierr); 31569ca1f37SDmitry Karpeev } 316a35b7b57SDmitry Karpeev for (d = 0; d < num_subdomains; ++d) { 317a35b7b57SDmitry Karpeev if (subdomain_names) {ierr = PetscFree(subdomain_names[d]);CHKERRQ(ierr);} 318a35b7b57SDmitry Karpeev if (inner_subdomain_is) {ierr = ISDestroy(&inner_subdomain_is[d]);CHKERRQ(ierr);} 319a35b7b57SDmitry Karpeev if (outer_subdomain_is) {ierr = ISDestroy(&outer_subdomain_is[d]);CHKERRQ(ierr);} 320fdc48646SDmitry Karpeev } 321a35b7b57SDmitry Karpeev ierr = PetscFree(subdomain_names);CHKERRQ(ierr); 322a35b7b57SDmitry Karpeev ierr = PetscFree(inner_subdomain_is);CHKERRQ(ierr); 323a35b7b57SDmitry Karpeev ierr = PetscFree(outer_subdomain_is);CHKERRQ(ierr); 324fdc48646SDmitry Karpeev } 32506b43e7eSDmitry Karpeev if (osm->n == PETSC_DECIDE) { /* still no subdomains; use one per processor */ 32644fe18e5SDmitry Karpeev osm->nmax = osm->n = 1; 327b1a0a8a3SJed Brown ierr = MPI_Comm_size(((PetscObject)pc)->comm,&size);CHKERRQ(ierr); 328f746d493SDmitry Karpeev osm->N = size; 329fdc48646SDmitry Karpeev } 33006b43e7eSDmitry Karpeev } 331a35b7b57SDmitry Karpeev if (!osm->iis) { 332a35b7b57SDmitry Karpeev /* 333a35b7b57SDmitry Karpeev The local number of subdomains was set just above, or in PCGASMSetTotalSubdomains(), or in PCGASMSetSubdomains(), 334a35b7b57SDmitry Karpeev but the actual subdomains have not been supplied (in PCGASMSetSubdomains()). 335a35b7b57SDmitry Karpeev We create the requisite number of inner subdomains on PETSC_COMM_SELF (for now). 336a35b7b57SDmitry Karpeev */ 337a35b7b57SDmitry Karpeev ierr = PCGASMCreateLocalSubdomains(pc->pmat,osm->overlap,osm->n,&osm->iis,&osm->ois);CHKERRQ(ierr); 338a35b7b57SDmitry Karpeev } 33906b43e7eSDmitry Karpeev if (osm->N == PETSC_DECIDE) { 3406ac3741eSJed Brown struct {PetscInt max,sum;} inwork,outwork; 341f746d493SDmitry Karpeev /* determine global number of subdomains and the max number of local subdomains */ 3426ac3741eSJed Brown inwork.max = osm->n; 3436ac3741eSJed Brown inwork.sum = osm->n; 3446ac3741eSJed Brown ierr = MPI_Allreduce(&inwork,&outwork,1,MPIU_2INT,PetscMaxSum_Op,((PetscObject)pc)->comm);CHKERRQ(ierr); 3456ac3741eSJed Brown osm->nmax = outwork.max; 3466ac3741eSJed Brown osm->N = outwork.sum; 347b1a0a8a3SJed Brown } 3486a4f0f73SDmitry Karpeev 349b1a0a8a3SJed Brown ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 350b1a0a8a3SJed Brown flg = PETSC_FALSE; 351f746d493SDmitry Karpeev ierr = PetscOptionsGetBool(prefix,"-pc_gasm_print_subdomains",&flg,PETSC_NULL);CHKERRQ(ierr); 352f746d493SDmitry Karpeev if (flg) { ierr = PCGASMPrintSubdomains(pc);CHKERRQ(ierr); } 353b1a0a8a3SJed Brown 354b1a0a8a3SJed Brown if (osm->sort_indices) { 355f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 3566a4f0f73SDmitry Karpeev ierr = ISSort(osm->ois[i]);CHKERRQ(ierr); 3576a4f0f73SDmitry Karpeev ierr = ISSort(osm->iis[i]);CHKERRQ(ierr); 358b1a0a8a3SJed Brown } 359b1a0a8a3SJed Brown } 3606a4f0f73SDmitry Karpeev /* 3616a4f0f73SDmitry Karpeev Merge the ISs, create merged vectors and restrictions. 3626a4f0f73SDmitry Karpeev */ 3636a4f0f73SDmitry Karpeev /* Merge outer subdomain ISs and construct a restriction onto the disjoint union of local outer subdomains. */ 3646a4f0f73SDmitry Karpeev on = 0; 365f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 3666a4f0f73SDmitry Karpeev ierr = ISGetLocalSize(osm->ois[i],&oni);CHKERRQ(ierr); 3676a4f0f73SDmitry Karpeev on += oni; 368f746d493SDmitry Karpeev } 3696a4f0f73SDmitry Karpeev ierr = PetscMalloc(on*sizeof(PetscInt), &oidx);CHKERRQ(ierr); 3706a4f0f73SDmitry Karpeev on = 0; 371f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 3726a4f0f73SDmitry Karpeev ierr = ISGetLocalSize(osm->ois[i],&oni);CHKERRQ(ierr); 3736a4f0f73SDmitry Karpeev ierr = ISGetIndices(osm->ois[i],&oidxi);CHKERRQ(ierr); 3746a4f0f73SDmitry Karpeev ierr = PetscMemcpy(oidx+on, oidxi, sizeof(PetscInt)*oni);CHKERRQ(ierr); 3756a4f0f73SDmitry Karpeev ierr = ISRestoreIndices(osm->ois[i], &oidxi);CHKERRQ(ierr); 3766a4f0f73SDmitry Karpeev on += oni; 377f746d493SDmitry Karpeev } 3786a4f0f73SDmitry Karpeev ierr = ISCreateGeneral(((PetscObject)(pc))->comm, on, oidx, PETSC_OWN_POINTER, &gois);CHKERRQ(ierr); 379f746d493SDmitry Karpeev ierr = MatGetVecs(pc->pmat,&x,&y);CHKERRQ(ierr); 3806a4f0f73SDmitry Karpeev ierr = VecCreateMPI(((PetscObject)pc)->comm, on, PETSC_DECIDE, &osm->gx);CHKERRQ(ierr); 381f746d493SDmitry Karpeev ierr = VecDuplicate(osm->gx,&osm->gy);CHKERRQ(ierr); 3826a4f0f73SDmitry Karpeev ierr = VecGetOwnershipRange(osm->gx, &gofirst, PETSC_NULL);CHKERRQ(ierr); 3836a4f0f73SDmitry Karpeev ierr = ISCreateStride(((PetscObject)pc)->comm,on,gofirst,1, &goid);CHKERRQ(ierr); 3846a4f0f73SDmitry Karpeev ierr = VecScatterCreate(x,gois,osm->gx,goid, &(osm->gorestriction));CHKERRQ(ierr); 3856a4f0f73SDmitry Karpeev ierr = VecDestroy(&x);CHKERRQ(ierr); 3867105d80bSDmitry Karpeev ierr = ISDestroy(&gois);CHKERRQ(ierr); 3876a4f0f73SDmitry Karpeev /* Merge inner subdomain ISs and construct a restriction onto the disjoint union of local inner subdomains. */ 3886a4f0f73SDmitry Karpeev { PetscInt ini; /* Number of indices the i-th a local inner subdomain. */ 3896a4f0f73SDmitry Karpeev PetscInt in; /* Number of indices in the disjoint uniont of local inner subdomains. */ 3906a4f0f73SDmitry Karpeev PetscInt *iidx; /* Global indices in the merged local inner subdomain. */ 3916a4f0f73SDmitry Karpeev PetscInt *ioidx; /* Global indices of the disjoint union of inner subdomains within the disjoint union of outer subdomains. */ 3926a4f0f73SDmitry Karpeev IS giis; /* IS for the disjoint union of inner subdomains. */ 3936a4f0f73SDmitry Karpeev IS giois; /* IS for the disjoint union of inner subdomains within the disjoint union of outer subdomains. */ 394f746d493SDmitry Karpeev /**/ 3956a4f0f73SDmitry Karpeev in = 0; 396f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 3976a4f0f73SDmitry Karpeev ierr = ISGetLocalSize(osm->iis[i],&ini);CHKERRQ(ierr); 3986a4f0f73SDmitry Karpeev in += ini; 399f746d493SDmitry Karpeev } 4006a4f0f73SDmitry Karpeev ierr = PetscMalloc(in*sizeof(PetscInt), &iidx);CHKERRQ(ierr); 4016a4f0f73SDmitry Karpeev ierr = PetscMalloc(in*sizeof(PetscInt), &ioidx);CHKERRQ(ierr); 4026a4f0f73SDmitry Karpeev ierr = VecGetOwnershipRange(osm->gx,&gofirst, PETSC_NULL);CHKERRQ(ierr); 4036a4f0f73SDmitry Karpeev in = 0; 4046a4f0f73SDmitry Karpeev on = 0; 405f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 4066a4f0f73SDmitry Karpeev const PetscInt *iidxi; /* Global indices of the i-th local inner subdomain. */ 4076a4f0f73SDmitry Karpeev ISLocalToGlobalMapping ltogi; /* Map from global to local indices of the i-th outer local subdomain. */ 4086a4f0f73SDmitry Karpeev PetscInt *ioidxi; /* Local indices of the i-th local inner subdomain within the local outer subdomain. */ 4096a4f0f73SDmitry Karpeev PetscInt ioni; /* Number of indices in ioidxi; if ioni != ini the inner subdomain is not a subdomain of the outer subdomain (error). */ 4106a4f0f73SDmitry Karpeev PetscInt k; 4116a4f0f73SDmitry Karpeev ierr = ISGetLocalSize(osm->iis[i],&ini);CHKERRQ(ierr); 4126a4f0f73SDmitry Karpeev ierr = ISGetLocalSize(osm->ois[i],&oni);CHKERRQ(ierr); 4136a4f0f73SDmitry Karpeev ierr = ISGetIndices(osm->iis[i],&iidxi);CHKERRQ(ierr); 4146a4f0f73SDmitry Karpeev ierr = PetscMemcpy(iidx+in, iidxi, sizeof(PetscInt)*ini);CHKERRQ(ierr); 4156a4f0f73SDmitry Karpeev ierr = ISLocalToGlobalMappingCreateIS(osm->ois[i],<ogi);CHKERRQ(ierr); 4166a4f0f73SDmitry Karpeev ioidxi = ioidx+in; 4176a4f0f73SDmitry Karpeev ierr = ISGlobalToLocalMappingApply(ltogi,IS_GTOLM_DROP,ini,iidxi,&ioni,ioidxi);CHKERRQ(ierr); 4186a4f0f73SDmitry Karpeev ierr = ISLocalToGlobalMappingDestroy(<ogi);CHKERRQ(ierr); 4196a4f0f73SDmitry Karpeev ierr = ISRestoreIndices(osm->iis[i], &iidxi);CHKERRQ(ierr); 4206a4f0f73SDmitry Karpeev if (ioni != ini) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_PLIB,"Inner subdomain %D contains %D indices outside of its outer subdomain", i, ini - ioni); 4216a4f0f73SDmitry Karpeev for (k = 0; k < ini; ++k) { 4226a4f0f73SDmitry Karpeev ioidxi[k] += gofirst+on; 423f746d493SDmitry Karpeev } 4246a4f0f73SDmitry Karpeev in += ini; 4256a4f0f73SDmitry Karpeev on += oni; 426f746d493SDmitry Karpeev } 4276a4f0f73SDmitry Karpeev ierr = ISCreateGeneral(((PetscObject)pc)->comm, in, iidx, PETSC_OWN_POINTER, &giis);CHKERRQ(ierr); 4286a4f0f73SDmitry Karpeev ierr = ISCreateGeneral(((PetscObject)pc)->comm, in, ioidx, PETSC_OWN_POINTER, &giois);CHKERRQ(ierr); 4296a4f0f73SDmitry Karpeev ierr = VecScatterCreate(y,giis,osm->gy,giois,&osm->girestriction);CHKERRQ(ierr); 4306a4f0f73SDmitry Karpeev ierr = VecDestroy(&y);CHKERRQ(ierr); 4316a4f0f73SDmitry Karpeev ierr = ISDestroy(&giis);CHKERRQ(ierr); 4326a4f0f73SDmitry Karpeev ierr = ISDestroy(&giois);CHKERRQ(ierr); 433b1a0a8a3SJed Brown } 4346a4f0f73SDmitry Karpeev ierr = ISDestroy(&goid);CHKERRQ(ierr); 435f746d493SDmitry Karpeev /* Create the subdomain work vectors. */ 436f746d493SDmitry Karpeev ierr = PetscMalloc(osm->n*sizeof(Vec),&osm->x);CHKERRQ(ierr); 437f746d493SDmitry Karpeev ierr = PetscMalloc(osm->n*sizeof(Vec),&osm->y);CHKERRQ(ierr); 438f746d493SDmitry Karpeev ierr = VecGetArray(osm->gx, &gxarray);CHKERRQ(ierr); 439f746d493SDmitry Karpeev ierr = VecGetArray(osm->gy, &gyarray);CHKERRQ(ierr); 4406a4f0f73SDmitry Karpeev for (i=0, on=0; i<osm->n; ++i, on += oni) { 4416a4f0f73SDmitry Karpeev PetscInt oNi; 4426a4f0f73SDmitry Karpeev ierr = ISGetLocalSize(osm->ois[i],&oni);CHKERRQ(ierr); 4436a4f0f73SDmitry Karpeev ierr = ISGetSize(osm->ois[i],&oNi);CHKERRQ(ierr); 4446a4f0f73SDmitry Karpeev ierr = VecCreateMPIWithArray(((PetscObject)(osm->ois[i]))->comm,1,oni,oNi,gxarray+on,&osm->x[i]);CHKERRQ(ierr); 4456a4f0f73SDmitry Karpeev ierr = VecCreateMPIWithArray(((PetscObject)(osm->ois[i]))->comm,1,oni,oNi,gyarray+on,&osm->y[i]);CHKERRQ(ierr); 446b1a0a8a3SJed Brown } 447f746d493SDmitry Karpeev ierr = VecRestoreArray(osm->gx, &gxarray);CHKERRQ(ierr); 448f746d493SDmitry Karpeev ierr = VecRestoreArray(osm->gy, &gyarray);CHKERRQ(ierr); 449b1a0a8a3SJed Brown /* Create the local solvers */ 450f746d493SDmitry Karpeev ierr = PetscMalloc(osm->n*sizeof(KSP *),&osm->ksp);CHKERRQ(ierr); 45144fe18e5SDmitry Karpeev for (i=0; i<osm->n; i++) { /* KSPs are local */ 4526a4f0f73SDmitry Karpeev ierr = KSPCreate(((PetscObject)(osm->ois[i]))->comm,&ksp);CHKERRQ(ierr); 453b1a0a8a3SJed Brown ierr = PetscLogObjectParent(pc,ksp);CHKERRQ(ierr); 454b1a0a8a3SJed Brown ierr = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,1);CHKERRQ(ierr); 455b1a0a8a3SJed Brown ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr); 456b1a0a8a3SJed Brown ierr = KSPGetPC(ksp,&subpc);CHKERRQ(ierr); 457b1a0a8a3SJed Brown ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 458b1a0a8a3SJed Brown ierr = KSPSetOptionsPrefix(ksp,prefix);CHKERRQ(ierr); 459b1a0a8a3SJed Brown ierr = KSPAppendOptionsPrefix(ksp,"sub_");CHKERRQ(ierr); 460b1a0a8a3SJed Brown osm->ksp[i] = ksp; 461b1a0a8a3SJed Brown } 462b1a0a8a3SJed Brown scall = MAT_INITIAL_MATRIX; 463b1a0a8a3SJed Brown 4646a4f0f73SDmitry Karpeev }/*if (!pc->setupcalled)*/ 4656a4f0f73SDmitry Karpeev else { 466b1a0a8a3SJed Brown /* 467b1a0a8a3SJed Brown Destroy the blocks from the previous iteration 468b1a0a8a3SJed Brown */ 469b1a0a8a3SJed Brown if (pc->flag == DIFFERENT_NONZERO_PATTERN) { 470f746d493SDmitry Karpeev ierr = MatDestroyMatrices(osm->n,&osm->pmat);CHKERRQ(ierr); 471b1a0a8a3SJed Brown scall = MAT_INITIAL_MATRIX; 472b1a0a8a3SJed Brown } 473b1a0a8a3SJed Brown } 474b1a0a8a3SJed Brown 475b1a0a8a3SJed Brown /* 476f746d493SDmitry Karpeev Extract out the submatrices. 477b1a0a8a3SJed Brown */ 47881a5c4aaSDmitry Karpeev if (size > 1) { 4796a4f0f73SDmitry Karpeev ierr = MatGetSubMatricesParallel(pc->pmat,osm->n,osm->ois, osm->ois,scall,&osm->pmat);CHKERRQ(ierr); 48081a5c4aaSDmitry Karpeev } 48181a5c4aaSDmitry Karpeev else { 4826a4f0f73SDmitry Karpeev ierr = MatGetSubMatrices(pc->pmat,osm->n,osm->ois, osm->ois,scall,&osm->pmat);CHKERRQ(ierr); 48381a5c4aaSDmitry Karpeev } 484b1a0a8a3SJed Brown if (scall == MAT_INITIAL_MATRIX) { 485b1a0a8a3SJed Brown ierr = PetscObjectGetOptionsPrefix((PetscObject)pc->pmat,&pprefix);CHKERRQ(ierr); 486f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 487b1a0a8a3SJed Brown ierr = PetscLogObjectParent(pc,osm->pmat[i]);CHKERRQ(ierr); 488b1a0a8a3SJed Brown ierr = PetscObjectSetOptionsPrefix((PetscObject)osm->pmat[i],pprefix);CHKERRQ(ierr); 489b1a0a8a3SJed Brown } 490b1a0a8a3SJed Brown } 491b1a0a8a3SJed Brown 492b1a0a8a3SJed Brown /* Return control to the user so that the submatrices can be modified (e.g., to apply 493b1a0a8a3SJed Brown different boundary conditions for the submatrices than for the global problem) */ 4946a4f0f73SDmitry Karpeev ierr = PCModifySubMatrices(pc,osm->n,osm->ois,osm->ois,osm->pmat,pc->modifysubmatricesP);CHKERRQ(ierr); 495b1a0a8a3SJed Brown 496b1a0a8a3SJed Brown /* 4976a4f0f73SDmitry Karpeev Loop over submatrices putting them into local ksps 498b1a0a8a3SJed Brown */ 499f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 500b1a0a8a3SJed Brown ierr = KSPSetOperators(osm->ksp[i],osm->pmat[i],osm->pmat[i],pc->flag);CHKERRQ(ierr); 501b1a0a8a3SJed Brown if (!pc->setupcalled) { 502b1a0a8a3SJed Brown ierr = KSPSetFromOptions(osm->ksp[i]);CHKERRQ(ierr); 503b1a0a8a3SJed Brown } 504b1a0a8a3SJed Brown } 505b1a0a8a3SJed Brown 506b1a0a8a3SJed Brown PetscFunctionReturn(0); 507b1a0a8a3SJed Brown } 508b1a0a8a3SJed Brown 509b1a0a8a3SJed Brown #undef __FUNCT__ 510f746d493SDmitry Karpeev #define __FUNCT__ "PCSetUpOnBlocks_GASM" 511f746d493SDmitry Karpeev static PetscErrorCode PCSetUpOnBlocks_GASM(PC pc) 512b1a0a8a3SJed Brown { 513f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 514b1a0a8a3SJed Brown PetscErrorCode ierr; 515b1a0a8a3SJed Brown PetscInt i; 516b1a0a8a3SJed Brown 517b1a0a8a3SJed Brown PetscFunctionBegin; 518f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 519b1a0a8a3SJed Brown ierr = KSPSetUp(osm->ksp[i]);CHKERRQ(ierr); 520b1a0a8a3SJed Brown } 521b1a0a8a3SJed Brown PetscFunctionReturn(0); 522b1a0a8a3SJed Brown } 523b1a0a8a3SJed Brown 524b1a0a8a3SJed Brown #undef __FUNCT__ 525f746d493SDmitry Karpeev #define __FUNCT__ "PCApply_GASM" 526f746d493SDmitry Karpeev static PetscErrorCode PCApply_GASM(PC pc,Vec x,Vec y) 527b1a0a8a3SJed Brown { 528f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 529b1a0a8a3SJed Brown PetscErrorCode ierr; 530f746d493SDmitry Karpeev PetscInt i; 531b1a0a8a3SJed Brown ScatterMode forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE; 532b1a0a8a3SJed Brown 533b1a0a8a3SJed Brown PetscFunctionBegin; 534b1a0a8a3SJed Brown /* 5356a4f0f73SDmitry Karpeev Support for limiting the restriction or interpolation only to the inner 536b1a0a8a3SJed Brown subdomain values (leaving the other values 0). 537b1a0a8a3SJed Brown */ 538f746d493SDmitry Karpeev if (!(osm->type & PC_GASM_RESTRICT)) { 539b1a0a8a3SJed Brown /* have to zero the work RHS since scatter may leave some slots empty */ 540f746d493SDmitry Karpeev ierr = VecZeroEntries(osm->gx);CHKERRQ(ierr); 5416a4f0f73SDmitry Karpeev ierr = VecScatterBegin(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 542b1a0a8a3SJed Brown } 5436a4f0f73SDmitry Karpeev else { 5446a4f0f73SDmitry Karpeev ierr = VecScatterBegin(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 545b1a0a8a3SJed Brown } 5466a4f0f73SDmitry Karpeev ierr = VecZeroEntries(osm->gy);CHKERRQ(ierr); 5476a4f0f73SDmitry Karpeev if (!(osm->type & PC_GASM_RESTRICT)) { 5486a4f0f73SDmitry Karpeev ierr = VecScatterEnd(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 5496a4f0f73SDmitry Karpeev } 5506a4f0f73SDmitry Karpeev else { 5516a4f0f73SDmitry Karpeev ierr = VecScatterEnd(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 5526a4f0f73SDmitry Karpeev } 55386b96d36SDmitry Karpeev /* do the subdomain solves */ 55486b96d36SDmitry Karpeev for (i=0; i<osm->n; ++i) { 555b1a0a8a3SJed Brown ierr = KSPSolve(osm->ksp[i],osm->x[i],osm->y[i]);CHKERRQ(ierr); 556b1a0a8a3SJed Brown } 5576a4f0f73SDmitry Karpeev ierr = VecZeroEntries(y);CHKERRQ(ierr); 5586a4f0f73SDmitry Karpeev if (!(osm->type & PC_GASM_INTERPOLATE)) { 5596a4f0f73SDmitry Karpeev ierr = VecScatterBegin(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); 5606a4f0f73SDmitry Karpeev ierr = VecScatterEnd(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); PetscFunctionReturn(0); 5616a4f0f73SDmitry Karpeev } 5626a4f0f73SDmitry Karpeev else { 5636a4f0f73SDmitry Karpeev ierr = VecScatterBegin(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); 5646a4f0f73SDmitry Karpeev ierr = VecScatterEnd(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); PetscFunctionReturn(0); 5656a4f0f73SDmitry Karpeev } 566b1a0a8a3SJed Brown } 567b1a0a8a3SJed Brown 568b1a0a8a3SJed Brown #undef __FUNCT__ 569f746d493SDmitry Karpeev #define __FUNCT__ "PCApplyTranspose_GASM" 570f746d493SDmitry Karpeev static PetscErrorCode PCApplyTranspose_GASM(PC pc,Vec x,Vec y) 571b1a0a8a3SJed Brown { 572f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 573b1a0a8a3SJed Brown PetscErrorCode ierr; 574f746d493SDmitry Karpeev PetscInt i; 575b1a0a8a3SJed Brown ScatterMode forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE; 576b1a0a8a3SJed Brown 577b1a0a8a3SJed Brown PetscFunctionBegin; 578b1a0a8a3SJed Brown /* 579b1a0a8a3SJed Brown Support for limiting the restriction or interpolation to only local 580b1a0a8a3SJed Brown subdomain values (leaving the other values 0). 581b1a0a8a3SJed Brown 582f746d493SDmitry Karpeev Note: these are reversed from the PCApply_GASM() because we are applying the 583b1a0a8a3SJed Brown transpose of the three terms 584b1a0a8a3SJed Brown */ 585f746d493SDmitry Karpeev if (!(osm->type & PC_GASM_INTERPOLATE)) { 586b1a0a8a3SJed Brown /* have to zero the work RHS since scatter may leave some slots empty */ 587f746d493SDmitry Karpeev ierr = VecZeroEntries(osm->gx);CHKERRQ(ierr); 5886a4f0f73SDmitry Karpeev ierr = VecScatterBegin(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 589b1a0a8a3SJed Brown } 5906a4f0f73SDmitry Karpeev else { 5916a4f0f73SDmitry Karpeev ierr = VecScatterBegin(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 592b1a0a8a3SJed Brown } 5936a4f0f73SDmitry Karpeev ierr = VecZeroEntries(osm->gy);CHKERRQ(ierr); 5946a4f0f73SDmitry Karpeev if (!(osm->type & PC_GASM_INTERPOLATE)) { 5956a4f0f73SDmitry Karpeev ierr = VecScatterEnd(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 5966a4f0f73SDmitry Karpeev } 5976a4f0f73SDmitry Karpeev else { 5986a4f0f73SDmitry Karpeev ierr = VecScatterEnd(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 5996a4f0f73SDmitry Karpeev } 600b1a0a8a3SJed Brown /* do the local solves */ 601ab3e923fSDmitry 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. */ 602b1a0a8a3SJed Brown ierr = KSPSolveTranspose(osm->ksp[i],osm->x[i],osm->y[i]);CHKERRQ(ierr); 603b1a0a8a3SJed Brown } 6046a4f0f73SDmitry Karpeev ierr = VecZeroEntries(y);CHKERRQ(ierr); 6056a4f0f73SDmitry Karpeev if (!(osm->type & PC_GASM_RESTRICT)) { 6066a4f0f73SDmitry Karpeev ierr = VecScatterBegin(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); 6076a4f0f73SDmitry Karpeev ierr = VecScatterEnd(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); 6086a4f0f73SDmitry Karpeev } 6096a4f0f73SDmitry Karpeev else { 6106a4f0f73SDmitry Karpeev ierr = VecScatterBegin(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); 6116a4f0f73SDmitry Karpeev ierr = VecScatterEnd(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); 6126a4f0f73SDmitry Karpeev } 6136a4f0f73SDmitry Karpeev 614b1a0a8a3SJed Brown PetscFunctionReturn(0); 615b1a0a8a3SJed Brown } 616b1a0a8a3SJed Brown 617b1a0a8a3SJed Brown #undef __FUNCT__ 618a06653b4SBarry Smith #define __FUNCT__ "PCReset_GASM" 619a06653b4SBarry Smith static PetscErrorCode PCReset_GASM(PC pc) 620b1a0a8a3SJed Brown { 621f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 622b1a0a8a3SJed Brown PetscErrorCode ierr; 623b1a0a8a3SJed Brown PetscInt i; 624b1a0a8a3SJed Brown 625b1a0a8a3SJed Brown PetscFunctionBegin; 626b1a0a8a3SJed Brown if (osm->ksp) { 627f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 628a06653b4SBarry Smith ierr = KSPReset(osm->ksp[i]);CHKERRQ(ierr); 629b1a0a8a3SJed Brown } 630b1a0a8a3SJed Brown } 631b1a0a8a3SJed Brown if (osm->pmat) { 632f746d493SDmitry Karpeev if (osm->n > 0) { 633ab3e923fSDmitry Karpeev ierr = MatDestroyMatrices(osm->n,&osm->pmat);CHKERRQ(ierr); 634b1a0a8a3SJed Brown } 635b1a0a8a3SJed Brown } 636d34fcf5fSBarry Smith if (osm->x) { 637f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 638fcfd50ebSBarry Smith ierr = VecDestroy(&osm->x[i]);CHKERRQ(ierr); 639fcfd50ebSBarry Smith ierr = VecDestroy(&osm->y[i]);CHKERRQ(ierr); 640b1a0a8a3SJed Brown } 641d34fcf5fSBarry Smith } 6426bf464f9SBarry Smith ierr = VecDestroy(&osm->gx);CHKERRQ(ierr); 6436bf464f9SBarry Smith ierr = VecDestroy(&osm->gy);CHKERRQ(ierr); 644ab3e923fSDmitry Karpeev 6456a4f0f73SDmitry Karpeev ierr = VecScatterDestroy(&osm->gorestriction);CHKERRQ(ierr); 6466a4f0f73SDmitry Karpeev ierr = VecScatterDestroy(&osm->girestriction);CHKERRQ(ierr); 6476a4f0f73SDmitry Karpeev ierr = PCGASMDestroySubdomains(osm->n,osm->ois,osm->iis);CHKERRQ(ierr); 6486a4f0f73SDmitry Karpeev osm->ois = 0; 6496a4f0f73SDmitry Karpeev osm->iis = 0; 650a06653b4SBarry Smith PetscFunctionReturn(0); 651a06653b4SBarry Smith } 652a06653b4SBarry Smith 653a06653b4SBarry Smith #undef __FUNCT__ 654a06653b4SBarry Smith #define __FUNCT__ "PCDestroy_GASM" 655a06653b4SBarry Smith static PetscErrorCode PCDestroy_GASM(PC pc) 656a06653b4SBarry Smith { 657a06653b4SBarry Smith PC_GASM *osm = (PC_GASM*)pc->data; 658a06653b4SBarry Smith PetscErrorCode ierr; 659a06653b4SBarry Smith PetscInt i; 660a06653b4SBarry Smith 661a06653b4SBarry Smith PetscFunctionBegin; 662a06653b4SBarry Smith ierr = PCReset_GASM(pc);CHKERRQ(ierr); 663a06653b4SBarry Smith if (osm->ksp) { 664a06653b4SBarry Smith for (i=0; i<osm->n; i++) { 6656bf464f9SBarry Smith ierr = KSPDestroy(&osm->ksp[i]);CHKERRQ(ierr); 666a06653b4SBarry Smith } 667a06653b4SBarry Smith ierr = PetscFree(osm->ksp);CHKERRQ(ierr); 668a06653b4SBarry Smith } 669a06653b4SBarry Smith ierr = PetscFree(osm->x);CHKERRQ(ierr); 670a06653b4SBarry Smith ierr = PetscFree(osm->y);CHKERRQ(ierr); 671c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 672b1a0a8a3SJed Brown PetscFunctionReturn(0); 673b1a0a8a3SJed Brown } 674b1a0a8a3SJed Brown 675b1a0a8a3SJed Brown #undef __FUNCT__ 676f746d493SDmitry Karpeev #define __FUNCT__ "PCSetFromOptions_GASM" 677*a6dfd86eSKarl Rupp static PetscErrorCode PCSetFromOptions_GASM(PC pc) 678*a6dfd86eSKarl Rupp { 679f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 680b1a0a8a3SJed Brown PetscErrorCode ierr; 681b1a0a8a3SJed Brown PetscInt blocks,ovl; 682b1a0a8a3SJed Brown PetscBool symset,flg; 683f746d493SDmitry Karpeev PCGASMType gasmtype; 684b1a0a8a3SJed Brown 685b1a0a8a3SJed Brown PetscFunctionBegin; 686b1a0a8a3SJed Brown /* set the type to symmetric if matrix is symmetric */ 687b1a0a8a3SJed Brown if (!osm->type_set && pc->pmat) { 688b1a0a8a3SJed Brown ierr = MatIsSymmetricKnown(pc->pmat,&symset,&flg);CHKERRQ(ierr); 689f746d493SDmitry Karpeev if (symset && flg) { osm->type = PC_GASM_BASIC; } 690b1a0a8a3SJed Brown } 69144fe18e5SDmitry Karpeev ierr = PetscOptionsHead("Generalized additive Schwarz options");CHKERRQ(ierr); 6926a4f0f73SDmitry Karpeev ierr = PetscOptionsInt("-pc_gasm_total_subdomains","Total number of subdomains across communicator","PCGASMSetTotalSubdomains",osm->n,&blocks,&flg);CHKERRQ(ierr); 6936a4f0f73SDmitry Karpeev osm->create_local = PETSC_TRUE; 6946a4f0f73SDmitry Karpeev ierr = PetscOptionsBool("-pc_gasm_subdomains_create_local","Whether to make autocreated subdomains local (true by default)","PCGASMSetTotalSubdomains",osm->create_local,&osm->create_local,&flg);CHKERRQ(ierr); 6956a4f0f73SDmitry Karpeev if (!osm->create_local) SETERRQ(((PetscObject)pc)->comm, PETSC_ERR_SUP, "No support for autocreation of nonlocal subdomains yet."); 6966a4f0f73SDmitry Karpeev 6976a4f0f73SDmitry Karpeev if (flg) {ierr = PCGASMSetTotalSubdomains(pc,blocks,osm->create_local);CHKERRQ(ierr); } 69806b43e7eSDmitry Karpeev ierr = PetscOptionsInt("-pc_gasm_overlap","Number of overlapping degrees of freedom","PCGASMSetOverlap",osm->overlap,&ovl,&flg);CHKERRQ(ierr); 699f746d493SDmitry Karpeev if (flg) {ierr = PCGASMSetOverlap(pc,ovl);CHKERRQ(ierr); } 700b1a0a8a3SJed Brown flg = PETSC_FALSE; 701f746d493SDmitry Karpeev ierr = PetscOptionsEnum("-pc_gasm_type","Type of restriction/extension","PCGASMSetType",PCGASMTypes,(PetscEnum)osm->type,(PetscEnum*)&gasmtype,&flg);CHKERRQ(ierr); 702f746d493SDmitry Karpeev if (flg) {ierr = PCGASMSetType(pc,gasmtype);CHKERRQ(ierr); } 703b1a0a8a3SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 704b1a0a8a3SJed Brown PetscFunctionReturn(0); 705b1a0a8a3SJed Brown } 706b1a0a8a3SJed Brown 707b1a0a8a3SJed Brown /*------------------------------------------------------------------------------------*/ 708b1a0a8a3SJed Brown 709b1a0a8a3SJed Brown EXTERN_C_BEGIN 710b1a0a8a3SJed Brown #undef __FUNCT__ 71106b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMSetSubdomains_GASM" 712a35b7b57SDmitry Karpeev PetscErrorCode PCGASMSetSubdomains_GASM(PC pc,PetscInt n,IS iis[],IS ois[]) 713b1a0a8a3SJed Brown { 714f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 715b1a0a8a3SJed Brown PetscErrorCode ierr; 716b1a0a8a3SJed Brown PetscInt i; 717b1a0a8a3SJed Brown 718b1a0a8a3SJed Brown PetscFunctionBegin; 71969ca1f37SDmitry Karpeev if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Each process must have 1 or more subdomains, n = %D",n); 720a35b7b57SDmitry Karpeev if (pc->setupcalled && (n != osm->n || iis || ois)) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetSubdomains() should be called before calling PCSetUp()."); 721b1a0a8a3SJed Brown 722b1a0a8a3SJed Brown if (!pc->setupcalled) { 72393c1ec32SDmitry Karpeev osm->n = n; 7246a4f0f73SDmitry Karpeev osm->ois = 0; 7256a4f0f73SDmitry Karpeev osm->iis = 0; 726a35b7b57SDmitry Karpeev if (ois) { 727a35b7b57SDmitry Karpeev for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)ois[i]);CHKERRQ(ierr);} 728b1a0a8a3SJed Brown } 729a35b7b57SDmitry Karpeev if (iis) { 730a35b7b57SDmitry Karpeev for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)iis[i]);CHKERRQ(ierr);} 731b1a0a8a3SJed Brown } 7326a4f0f73SDmitry Karpeev ierr = PCGASMDestroySubdomains(osm->n,osm->iis,osm->ois);CHKERRQ(ierr); 733a35b7b57SDmitry Karpeev if (ois) { 7346a4f0f73SDmitry Karpeev ierr = PetscMalloc(n*sizeof(IS),&osm->ois);CHKERRQ(ierr); 735a35b7b57SDmitry Karpeev for (i=0; i<n; i++) { osm->ois[i] = ois[i]; } 7366a4f0f73SDmitry Karpeev /* Flag indicating that the user has set outer subdomains, so PCGASM should not increase their size. */ 737b1a0a8a3SJed Brown osm->overlap = -1; 738a35b7b57SDmitry Karpeev if (!iis) { 7396a4f0f73SDmitry Karpeev ierr = PetscMalloc(n*sizeof(IS),&osm->iis);CHKERRQ(ierr); 740a35b7b57SDmitry Karpeev for (i=0; i<n; i++) { 741a35b7b57SDmitry Karpeev for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)ois[i]);CHKERRQ(ierr);} 742a35b7b57SDmitry Karpeev osm->iis[i] = ois[i]; 743a35b7b57SDmitry Karpeev } 744a35b7b57SDmitry Karpeev } 745a35b7b57SDmitry Karpeev } 746a35b7b57SDmitry Karpeev if (iis) { 747a35b7b57SDmitry Karpeev ierr = PetscMalloc(n*sizeof(IS),&osm->iis);CHKERRQ(ierr); 748a35b7b57SDmitry Karpeev for (i=0; i<n; i++) { osm->iis[i] = iis[i]; } 749a35b7b57SDmitry Karpeev if (!ois) { 750a35b7b57SDmitry Karpeev ierr = PetscMalloc(n*sizeof(IS),&osm->ois);CHKERRQ(ierr); 751a35b7b57SDmitry Karpeev for (i=0; i<n; i++) { 752a35b7b57SDmitry Karpeev for (i=0; i<n; i++) { 753a35b7b57SDmitry Karpeev ierr = PetscObjectReference((PetscObject)iis[i]);CHKERRQ(ierr); 754a35b7b57SDmitry Karpeev osm->ois[i] = iis[i]; 755a35b7b57SDmitry Karpeev } 756a35b7b57SDmitry Karpeev } 757a35b7b57SDmitry Karpeev if (osm->overlap > 0) { 758a35b7b57SDmitry Karpeev /* Extend the "overlapping" regions by a number of steps */ 759a35b7b57SDmitry Karpeev ierr = MatIncreaseOverlap(pc->pmat,osm->n,osm->ois,osm->overlap);CHKERRQ(ierr); 760a35b7b57SDmitry Karpeev } 761a35b7b57SDmitry Karpeev } 762b1a0a8a3SJed Brown } 763b1a0a8a3SJed Brown } 764b1a0a8a3SJed Brown PetscFunctionReturn(0); 765b1a0a8a3SJed Brown } 766b1a0a8a3SJed Brown EXTERN_C_END 767b1a0a8a3SJed Brown 768b1a0a8a3SJed Brown EXTERN_C_BEGIN 769b1a0a8a3SJed Brown #undef __FUNCT__ 7706a4f0f73SDmitry Karpeev #define __FUNCT__ "PCGASMSetTotalSubdomains_GASM" 7710adebc6cSBarry Smith PetscErrorCode PCGASMSetTotalSubdomains_GASM(PC pc,PetscInt N, PetscBool create_local) 7720adebc6cSBarry Smith { 773f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 774b1a0a8a3SJed Brown PetscErrorCode ierr; 775b1a0a8a3SJed Brown PetscMPIInt rank,size; 776b1a0a8a3SJed Brown PetscInt n; 7776a4f0f73SDmitry Karpeev PetscInt Nmin, Nmax; 778b1a0a8a3SJed Brown PetscFunctionBegin; 7796a4f0f73SDmitry Karpeev if (!create_local) SETERRQ(((PetscObject)pc)->comm, PETSC_ERR_SUP, "No suppor for autocreation of nonlocal subdomains."); 7806a4f0f73SDmitry Karpeev if (N < 1) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Total number of subdomains must be > 0, N = %D",N); 781adcec1e5SJed Brown ierr = MPI_Allreduce(&N,&Nmin,1,MPIU_INT,MPIU_MIN,((PetscObject)pc)->comm);CHKERRQ(ierr); 782adcec1e5SJed Brown ierr = MPI_Allreduce(&N,&Nmax,1,MPIU_INT,MPIU_MAX,((PetscObject)pc)->comm);CHKERRQ(ierr); 7836a4f0f73SDmitry Karpeev if (Nmin != Nmax) 7846a4f0f73SDmitry Karpeev SETERRQ2(((PetscObject)pc)->comm, PETSC_ERR_ARG_WRONG, "All processors must use the same number of subdomains. min(N) = %D != %D = max(N)", Nmin, Nmax); 785b1a0a8a3SJed Brown 7866a4f0f73SDmitry Karpeev osm->create_local = create_local; 787b1a0a8a3SJed Brown /* 788b1a0a8a3SJed Brown Split the subdomains equally among all processors 789b1a0a8a3SJed Brown */ 790b1a0a8a3SJed Brown ierr = MPI_Comm_rank(((PetscObject)pc)->comm,&rank);CHKERRQ(ierr); 791b1a0a8a3SJed Brown ierr = MPI_Comm_size(((PetscObject)pc)->comm,&size);CHKERRQ(ierr); 792b1a0a8a3SJed Brown n = N/size + ((N % size) > rank); 7936a4f0f73SDmitry Karpeev if (!n) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Process %d must have at least one subdomain: total processors %d total blocks %D",(int)rank,(int)size,N); 794f746d493SDmitry Karpeev if (pc->setupcalled && n != osm->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetTotalSubdomains() should be called before PCSetUp()."); 795b1a0a8a3SJed Brown if (!pc->setupcalled) { 7966a4f0f73SDmitry Karpeev ierr = PCGASMDestroySubdomains(osm->n,osm->iis,osm->ois);CHKERRQ(ierr); 79793c1ec32SDmitry Karpeev osm->N = N; 798f746d493SDmitry Karpeev osm->n = n; 79906b43e7eSDmitry Karpeev osm->nmax = N/size + ((N%size)?1:0); 8006a4f0f73SDmitry Karpeev osm->ois = 0; 8016a4f0f73SDmitry Karpeev osm->iis = 0; 802b1a0a8a3SJed Brown } 803b1a0a8a3SJed Brown PetscFunctionReturn(0); 80406b43e7eSDmitry Karpeev } 805b1a0a8a3SJed Brown EXTERN_C_END 806b1a0a8a3SJed Brown 807b1a0a8a3SJed Brown EXTERN_C_BEGIN 808b1a0a8a3SJed Brown #undef __FUNCT__ 809f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetOverlap_GASM" 8107087cfbeSBarry Smith PetscErrorCode PCGASMSetOverlap_GASM(PC pc,PetscInt ovl) 811b1a0a8a3SJed Brown { 812f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 813b1a0a8a3SJed Brown 814b1a0a8a3SJed Brown PetscFunctionBegin; 815b1a0a8a3SJed Brown if (ovl < 0) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap value requested"); 816f746d493SDmitry Karpeev if (pc->setupcalled && ovl != osm->overlap) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetOverlap() should be called before PCSetUp()."); 817b1a0a8a3SJed Brown if (!pc->setupcalled) { 818b1a0a8a3SJed Brown osm->overlap = ovl; 819b1a0a8a3SJed Brown } 820b1a0a8a3SJed Brown PetscFunctionReturn(0); 821b1a0a8a3SJed Brown } 822b1a0a8a3SJed Brown EXTERN_C_END 823b1a0a8a3SJed Brown 824b1a0a8a3SJed Brown EXTERN_C_BEGIN 825b1a0a8a3SJed Brown #undef __FUNCT__ 826f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetType_GASM" 8277087cfbeSBarry Smith PetscErrorCode PCGASMSetType_GASM(PC pc,PCGASMType type) 828b1a0a8a3SJed Brown { 829f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 830b1a0a8a3SJed Brown 831b1a0a8a3SJed Brown PetscFunctionBegin; 832b1a0a8a3SJed Brown osm->type = type; 833b1a0a8a3SJed Brown osm->type_set = PETSC_TRUE; 834b1a0a8a3SJed Brown PetscFunctionReturn(0); 835b1a0a8a3SJed Brown } 836b1a0a8a3SJed Brown EXTERN_C_END 837b1a0a8a3SJed Brown 838b1a0a8a3SJed Brown EXTERN_C_BEGIN 839b1a0a8a3SJed Brown #undef __FUNCT__ 840f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetSortIndices_GASM" 8417087cfbeSBarry Smith PetscErrorCode PCGASMSetSortIndices_GASM(PC pc,PetscBool doSort) 842b1a0a8a3SJed Brown { 843f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 844b1a0a8a3SJed Brown 845b1a0a8a3SJed Brown PetscFunctionBegin; 846b1a0a8a3SJed Brown osm->sort_indices = doSort; 847b1a0a8a3SJed Brown PetscFunctionReturn(0); 848b1a0a8a3SJed Brown } 849b1a0a8a3SJed Brown EXTERN_C_END 850b1a0a8a3SJed Brown 851b1a0a8a3SJed Brown EXTERN_C_BEGIN 852b1a0a8a3SJed Brown #undef __FUNCT__ 853f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMGetSubKSP_GASM" 85444fe18e5SDmitry Karpeev /* 85544fe18e5SDmitry Karpeev FIX: This routine might need to be modified once multiple ranks per subdomain are allowed. 85644fe18e5SDmitry Karpeev In particular, it would upset the global subdomain number calculation. 85744fe18e5SDmitry Karpeev */ 8587087cfbeSBarry Smith PetscErrorCode PCGASMGetSubKSP_GASM(PC pc,PetscInt *n,PetscInt *first,KSP **ksp) 859b1a0a8a3SJed Brown { 860f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 861b1a0a8a3SJed Brown PetscErrorCode ierr; 862b1a0a8a3SJed Brown 863b1a0a8a3SJed Brown PetscFunctionBegin; 864ab3e923fSDmitry Karpeev if (osm->n < 1) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ORDER,"Need to call PCSetUP() on PC (or KSPSetUp() on the outer KSP object) before calling here"); 865b1a0a8a3SJed Brown 866ab3e923fSDmitry Karpeev if (n) { 867ab3e923fSDmitry Karpeev *n = osm->n; 868b1a0a8a3SJed Brown } 869ab3e923fSDmitry Karpeev if (first) { 870ab3e923fSDmitry Karpeev ierr = MPI_Scan(&osm->n,first,1,MPIU_INT,MPI_SUM,((PetscObject)pc)->comm);CHKERRQ(ierr); 871ab3e923fSDmitry Karpeev *first -= osm->n; 872b1a0a8a3SJed Brown } 873b1a0a8a3SJed Brown if (ksp) { 874b1a0a8a3SJed Brown /* Assume that local solves are now different; not necessarily 87506b43e7eSDmitry Karpeev true, though! This flag is used only for PCView_GASM() */ 876b1a0a8a3SJed Brown *ksp = osm->ksp; 8776a4f0f73SDmitry Karpeev osm->same_subdomain_solvers = PETSC_FALSE; 878b1a0a8a3SJed Brown } 879b1a0a8a3SJed Brown PetscFunctionReturn(0); 880ab3e923fSDmitry Karpeev }/* PCGASMGetSubKSP_GASM() */ 881b1a0a8a3SJed Brown EXTERN_C_END 882b1a0a8a3SJed Brown 883b1a0a8a3SJed Brown 884b1a0a8a3SJed Brown #undef __FUNCT__ 88506b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMSetSubdomains" 886b1a0a8a3SJed Brown /*@C 88706b43e7eSDmitry Karpeev PCGASMSetSubdomains - Sets the subdomains for this processor 88806b43e7eSDmitry Karpeev for the additive Schwarz preconditioner. 889b1a0a8a3SJed Brown 890b1a0a8a3SJed Brown Collective on PC 891b1a0a8a3SJed Brown 892b1a0a8a3SJed Brown Input Parameters: 893b1a0a8a3SJed Brown + pc - the preconditioner context 89406b43e7eSDmitry Karpeev . n - the number of subdomains for this processor 8956a4f0f73SDmitry Karpeev . iis - the index sets that define this processor's local inner subdomains 896b1a0a8a3SJed Brown (or PETSC_NULL for PETSc to determine subdomains) 8976a4f0f73SDmitry Karpeev - ois- the index sets that define this processor's local outer subdomains 8986a4f0f73SDmitry Karpeev (or PETSC_NULL to use the same as iis) 899b1a0a8a3SJed Brown 900b1a0a8a3SJed Brown Notes: 90106b43e7eSDmitry Karpeev The IS indices use the parallel, global numbering of the vector entries. 9026a4f0f73SDmitry Karpeev Inner subdomains are those where the correction is applied. 9036a4f0f73SDmitry Karpeev Outer subdomains are those where the residual necessary to obtain the 9046a4f0f73SDmitry Karpeev corrections is obtained (see PCGASMType for the use of inner/outer subdomains). 9056a4f0f73SDmitry Karpeev Both inner and outer subdomains can extend over several processors. 9066a4f0f73SDmitry Karpeev This processor's portion of a subdomain is known as a local subdomain. 9076a4f0f73SDmitry Karpeev 9086a4f0f73SDmitry Karpeev By default the GASM preconditioner uses 1 (local) subdomain per processor. 9096a4f0f73SDmitry Karpeev Use PCGASMSetTotalSubdomains() to set the total number of subdomains across 9106a4f0f73SDmitry Karpeev all processors that PCGASM will create automatically, and to specify whether 9116a4f0f73SDmitry Karpeev they should be local or not. 9126a4f0f73SDmitry Karpeev 913b1a0a8a3SJed Brown 914b1a0a8a3SJed Brown Level: advanced 915b1a0a8a3SJed Brown 91606b43e7eSDmitry Karpeev .keywords: PC, GASM, set, subdomains, additive Schwarz 917b1a0a8a3SJed Brown 9186a4f0f73SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetOverlap(), PCGASMGetSubKSP(), 91906b43e7eSDmitry Karpeev PCGASMCreateSubdomains2D(), PCGASMGetSubdomains() 920b1a0a8a3SJed Brown @*/ 9216a4f0f73SDmitry Karpeev PetscErrorCode PCGASMSetSubdomains(PC pc,PetscInt n,IS iis[],IS ois[]) 922b1a0a8a3SJed Brown { 923b1a0a8a3SJed Brown PetscErrorCode ierr; 924b1a0a8a3SJed Brown 925b1a0a8a3SJed Brown PetscFunctionBegin; 926b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 9276a4f0f73SDmitry Karpeev ierr = PetscTryMethod(pc,"PCGASMSetSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,n,iis,ois));CHKERRQ(ierr); 928b1a0a8a3SJed Brown PetscFunctionReturn(0); 929b1a0a8a3SJed Brown } 930b1a0a8a3SJed Brown 931b1a0a8a3SJed Brown #undef __FUNCT__ 9326a4f0f73SDmitry Karpeev #define __FUNCT__ "PCGASMSetTotalSubdomains" 933b1a0a8a3SJed Brown /*@C 9346a4f0f73SDmitry Karpeev PCGASMSetTotalSubdomains - Sets the total number of subdomains to use in the generalized additive 9356a4f0f73SDmitry Karpeev Schwarz preconditioner. The number of subdomains is cumulative across all processors in pc's 9366a4f0f73SDmitry Karpeev communicator. Either all or no processors in the PC communicator must call this routine with 9376a4f0f73SDmitry Karpeev the same N. The subdomains will be created automatically during PCSetUp(). 938b1a0a8a3SJed Brown 939b1a0a8a3SJed Brown Collective on PC 940b1a0a8a3SJed Brown 941b1a0a8a3SJed Brown Input Parameters: 942b1a0a8a3SJed Brown + pc - the preconditioner context 9436a4f0f73SDmitry Karpeev . N - the total number of subdomains cumulative across all processors 9446a4f0f73SDmitry Karpeev - create_local - whether the subdomains to be created are to be local 945b1a0a8a3SJed Brown 946b1a0a8a3SJed Brown Options Database Key: 9476a4f0f73SDmitry Karpeev To set the total number of subdomains and let PCGASM autocreate them, rather than specify the index sets, use the following options: 9486a4f0f73SDmitry Karpeev + -pc_gasm_total_subdomains <n> - sets the total number of subdomains to be autocreated by PCGASM 9496a4f0f73SDmitry Karpeev - -pc_gasm_subdomains_create_local <true|false> - whether autocreated subdomains should be local or not (default is true) 950b1a0a8a3SJed Brown 95106b43e7eSDmitry Karpeev By default the GASM preconditioner uses 1 subdomain per processor. 952b1a0a8a3SJed Brown 9536a4f0f73SDmitry Karpeev 9546a4f0f73SDmitry Karpeev Use PCGASMSetSubdomains() to set subdomains explicitly or to set different numbers 9556a4f0f73SDmitry Karpeev of subdomains per processor. 956b1a0a8a3SJed Brown 957b1a0a8a3SJed Brown Level: advanced 958b1a0a8a3SJed Brown 959f746d493SDmitry Karpeev .keywords: PC, GASM, set, total, global, subdomains, additive Schwarz 960b1a0a8a3SJed Brown 96106b43e7eSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMSetOverlap(), PCGASMGetSubKSP(), 962f746d493SDmitry Karpeev PCGASMCreateSubdomains2D() 963b1a0a8a3SJed Brown @*/ 9646a4f0f73SDmitry Karpeev PetscErrorCode PCGASMSetTotalSubdomains(PC pc,PetscInt N, PetscBool create_local) 965b1a0a8a3SJed Brown { 966b1a0a8a3SJed Brown PetscErrorCode ierr; 967b1a0a8a3SJed Brown 968b1a0a8a3SJed Brown PetscFunctionBegin; 969b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 9706a4f0f73SDmitry Karpeev ierr = PetscTryMethod(pc,"PCGASMSetTotalSubdomains_C",(PC,PetscInt,PetscBool),(pc,N,create_local));CHKERRQ(ierr); 971b1a0a8a3SJed Brown PetscFunctionReturn(0); 972b1a0a8a3SJed Brown } 973b1a0a8a3SJed Brown 974b1a0a8a3SJed Brown #undef __FUNCT__ 975f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetOverlap" 976b1a0a8a3SJed Brown /*@ 977f746d493SDmitry Karpeev PCGASMSetOverlap - Sets the overlap between a pair of subdomains for the 978b1a0a8a3SJed Brown additive Schwarz preconditioner. Either all or no processors in the 979b1a0a8a3SJed Brown PC communicator must call this routine. 980b1a0a8a3SJed Brown 981b1a0a8a3SJed Brown Logically Collective on PC 982b1a0a8a3SJed Brown 983b1a0a8a3SJed Brown Input Parameters: 984b1a0a8a3SJed Brown + pc - the preconditioner context 985b1a0a8a3SJed Brown - ovl - the amount of overlap between subdomains (ovl >= 0, default value = 1) 986b1a0a8a3SJed Brown 987b1a0a8a3SJed Brown Options Database Key: 98806b43e7eSDmitry Karpeev . -pc_gasm_overlap <overlap> - Sets overlap 989b1a0a8a3SJed Brown 990b1a0a8a3SJed Brown Notes: 99106b43e7eSDmitry Karpeev By default the GASM preconditioner uses 1 subdomain per processor. To use 9926a4f0f73SDmitry Karpeev multiple subdomain per perocessor, see PCGASMSetTotalSubdomains() or 9936a4f0f73SDmitry Karpeev PCGASMSetSubdomains() (and the option -pc_gasm_total_subdomains <n>). 994b1a0a8a3SJed Brown 995b1a0a8a3SJed Brown The overlap defaults to 1, so if one desires that no additional 996b1a0a8a3SJed Brown overlap be computed beyond what may have been set with a call to 9976a4f0f73SDmitry Karpeev PCGASMSetTotalSubdomains() or PCGASMSetSubdomains(), then ovl 998b1a0a8a3SJed Brown must be set to be 0. In particular, if one does not explicitly set 99906b43e7eSDmitry Karpeev the subdomains in application code, then all overlap would be computed 1000f746d493SDmitry Karpeev internally by PETSc, and using an overlap of 0 would result in an GASM 1001b1a0a8a3SJed Brown variant that is equivalent to the block Jacobi preconditioner. 1002b1a0a8a3SJed Brown 1003b1a0a8a3SJed Brown Note that one can define initial index sets with any overlap via 100406b43e7eSDmitry Karpeev PCGASMSetSubdomains(); the routine PCGASMSetOverlap() merely allows 100506b43e7eSDmitry Karpeev PETSc to extend that overlap further, if desired. 1006b1a0a8a3SJed Brown 1007b1a0a8a3SJed Brown Level: intermediate 1008b1a0a8a3SJed Brown 1009f746d493SDmitry Karpeev .keywords: PC, GASM, set, overlap 1010b1a0a8a3SJed Brown 10116a4f0f73SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetSubdomains(), PCGASMGetSubKSP(), 101206b43e7eSDmitry Karpeev PCGASMCreateSubdomains2D(), PCGASMGetSubdomains() 1013b1a0a8a3SJed Brown @*/ 10147087cfbeSBarry Smith PetscErrorCode PCGASMSetOverlap(PC pc,PetscInt ovl) 1015b1a0a8a3SJed Brown { 1016b1a0a8a3SJed Brown PetscErrorCode ierr; 1017b1a0a8a3SJed Brown 1018b1a0a8a3SJed Brown PetscFunctionBegin; 1019b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1020b1a0a8a3SJed Brown PetscValidLogicalCollectiveInt(pc,ovl,2); 1021f746d493SDmitry Karpeev ierr = PetscTryMethod(pc,"PCGASMSetOverlap_C",(PC,PetscInt),(pc,ovl));CHKERRQ(ierr); 1022b1a0a8a3SJed Brown PetscFunctionReturn(0); 1023b1a0a8a3SJed Brown } 1024b1a0a8a3SJed Brown 1025b1a0a8a3SJed Brown #undef __FUNCT__ 1026f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetType" 1027b1a0a8a3SJed Brown /*@ 1028f746d493SDmitry Karpeev PCGASMSetType - Sets the type of restriction and interpolation used 1029b1a0a8a3SJed Brown for local problems in the additive Schwarz method. 1030b1a0a8a3SJed Brown 1031b1a0a8a3SJed Brown Logically Collective on PC 1032b1a0a8a3SJed Brown 1033b1a0a8a3SJed Brown Input Parameters: 1034b1a0a8a3SJed Brown + pc - the preconditioner context 1035f746d493SDmitry Karpeev - type - variant of GASM, one of 1036b1a0a8a3SJed Brown .vb 1037f746d493SDmitry Karpeev PC_GASM_BASIC - full interpolation and restriction 1038f746d493SDmitry Karpeev PC_GASM_RESTRICT - full restriction, local processor interpolation 1039f746d493SDmitry Karpeev PC_GASM_INTERPOLATE - full interpolation, local processor restriction 1040f746d493SDmitry Karpeev PC_GASM_NONE - local processor restriction and interpolation 1041b1a0a8a3SJed Brown .ve 1042b1a0a8a3SJed Brown 1043b1a0a8a3SJed Brown Options Database Key: 1044f746d493SDmitry Karpeev . -pc_gasm_type [basic,restrict,interpolate,none] - Sets GASM type 1045b1a0a8a3SJed Brown 1046b1a0a8a3SJed Brown Level: intermediate 1047b1a0a8a3SJed Brown 1048f746d493SDmitry Karpeev .keywords: PC, GASM, set, type 1049b1a0a8a3SJed Brown 10506a4f0f73SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetSubdomains(), PCGASMGetSubKSP(), 1051f746d493SDmitry Karpeev PCGASMCreateSubdomains2D() 1052b1a0a8a3SJed Brown @*/ 10537087cfbeSBarry Smith PetscErrorCode PCGASMSetType(PC pc,PCGASMType type) 1054b1a0a8a3SJed Brown { 1055b1a0a8a3SJed Brown PetscErrorCode ierr; 1056b1a0a8a3SJed Brown 1057b1a0a8a3SJed Brown PetscFunctionBegin; 1058b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1059b1a0a8a3SJed Brown PetscValidLogicalCollectiveEnum(pc,type,2); 1060f746d493SDmitry Karpeev ierr = PetscTryMethod(pc,"PCGASMSetType_C",(PC,PCGASMType),(pc,type));CHKERRQ(ierr); 1061b1a0a8a3SJed Brown PetscFunctionReturn(0); 1062b1a0a8a3SJed Brown } 1063b1a0a8a3SJed Brown 1064b1a0a8a3SJed Brown #undef __FUNCT__ 1065f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetSortIndices" 1066b1a0a8a3SJed Brown /*@ 1067f746d493SDmitry Karpeev PCGASMSetSortIndices - Determines whether subdomain indices are sorted. 1068b1a0a8a3SJed Brown 1069b1a0a8a3SJed Brown Logically Collective on PC 1070b1a0a8a3SJed Brown 1071b1a0a8a3SJed Brown Input Parameters: 1072b1a0a8a3SJed Brown + pc - the preconditioner context 1073b1a0a8a3SJed Brown - doSort - sort the subdomain indices 1074b1a0a8a3SJed Brown 1075b1a0a8a3SJed Brown Level: intermediate 1076b1a0a8a3SJed Brown 1077f746d493SDmitry Karpeev .keywords: PC, GASM, set, type 1078b1a0a8a3SJed Brown 10796a4f0f73SDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMSetTotalSubdomains(), PCGASMGetSubKSP(), 1080f746d493SDmitry Karpeev PCGASMCreateSubdomains2D() 1081b1a0a8a3SJed Brown @*/ 10827087cfbeSBarry Smith PetscErrorCode PCGASMSetSortIndices(PC pc,PetscBool doSort) 1083b1a0a8a3SJed Brown { 1084b1a0a8a3SJed Brown PetscErrorCode ierr; 1085b1a0a8a3SJed Brown 1086b1a0a8a3SJed Brown PetscFunctionBegin; 1087b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1088b1a0a8a3SJed Brown PetscValidLogicalCollectiveBool(pc,doSort,2); 1089f746d493SDmitry Karpeev ierr = PetscTryMethod(pc,"PCGASMSetSortIndices_C",(PC,PetscBool),(pc,doSort));CHKERRQ(ierr); 1090b1a0a8a3SJed Brown PetscFunctionReturn(0); 1091b1a0a8a3SJed Brown } 1092b1a0a8a3SJed Brown 1093b1a0a8a3SJed Brown #undef __FUNCT__ 1094f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMGetSubKSP" 1095b1a0a8a3SJed Brown /*@C 1096f746d493SDmitry Karpeev PCGASMGetSubKSP - Gets the local KSP contexts for all blocks on 1097b1a0a8a3SJed Brown this processor. 1098b1a0a8a3SJed Brown 1099b1a0a8a3SJed Brown Collective on PC iff first_local is requested 1100b1a0a8a3SJed Brown 1101b1a0a8a3SJed Brown Input Parameter: 1102b1a0a8a3SJed Brown . pc - the preconditioner context 1103b1a0a8a3SJed Brown 1104b1a0a8a3SJed Brown Output Parameters: 1105b1a0a8a3SJed Brown + n_local - the number of blocks on this processor or PETSC_NULL 1106b1a0a8a3SJed Brown . first_local - the global number of the first block on this processor or PETSC_NULL, 1107b1a0a8a3SJed Brown all processors must request or all must pass PETSC_NULL 1108b1a0a8a3SJed Brown - ksp - the array of KSP contexts 1109b1a0a8a3SJed Brown 1110b1a0a8a3SJed Brown Note: 1111f746d493SDmitry Karpeev After PCGASMGetSubKSP() the array of KSPes is not to be freed 1112b1a0a8a3SJed Brown 1113b1a0a8a3SJed Brown Currently for some matrix implementations only 1 block per processor 1114b1a0a8a3SJed Brown is supported. 1115b1a0a8a3SJed Brown 1116f746d493SDmitry Karpeev You must call KSPSetUp() before calling PCGASMGetSubKSP(). 1117b1a0a8a3SJed Brown 1118b1a0a8a3SJed Brown Level: advanced 1119b1a0a8a3SJed Brown 1120f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, get, sub, KSP, context 1121b1a0a8a3SJed Brown 11226a4f0f73SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetSubdomains(), PCGASMSetOverlap(), 1123f746d493SDmitry Karpeev PCGASMCreateSubdomains2D(), 1124b1a0a8a3SJed Brown @*/ 11257087cfbeSBarry Smith PetscErrorCode PCGASMGetSubKSP(PC pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp[]) 1126b1a0a8a3SJed Brown { 1127b1a0a8a3SJed Brown PetscErrorCode ierr; 1128b1a0a8a3SJed Brown 1129b1a0a8a3SJed Brown PetscFunctionBegin; 1130b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1131f746d493SDmitry Karpeev ierr = PetscUseMethod(pc,"PCGASMGetSubKSP_C",(PC,PetscInt*,PetscInt*,KSP **),(pc,n_local,first_local,ksp));CHKERRQ(ierr); 1132b1a0a8a3SJed Brown PetscFunctionReturn(0); 1133b1a0a8a3SJed Brown } 1134b1a0a8a3SJed Brown 1135b1a0a8a3SJed Brown /* -------------------------------------------------------------------------------------*/ 1136b1a0a8a3SJed Brown /*MC 1137f746d493SDmitry Karpeev PCGASM - Use the (restricted) additive Schwarz method, each block is (approximately) solved with 1138b1a0a8a3SJed Brown its own KSP object. 1139b1a0a8a3SJed Brown 1140b1a0a8a3SJed Brown Options Database Keys: 114106b43e7eSDmitry Karpeev + -pc_gasm_total_block_count <n> - Sets total number of local subdomains (known as blocks) to be distributed among processors 114206b43e7eSDmitry Karpeev . -pc_gasm_view_subdomains - activates the printing of subdomain indices in PCView(), -ksp_view or -snes_view 114306b43e7eSDmitry Karpeev . -pc_gasm_print_subdomains - activates the printing of subdomain indices in PCSetUp() 114406b43e7eSDmitry Karpeev . -pc_gasm_overlap <ovl> - Sets overlap by which to (automatically) extend local subdomains 1145f746d493SDmitry Karpeev - -pc_gasm_type [basic,restrict,interpolate,none] - Sets GASM type 1146b1a0a8a3SJed Brown 1147b1a0a8a3SJed Brown IMPORTANT: If you run with, for example, 3 blocks on 1 processor or 3 blocks on 3 processors you 1148f746d493SDmitry Karpeev will get a different convergence rate due to the default option of -pc_gasm_type restrict. Use 1149f746d493SDmitry Karpeev -pc_gasm_type basic to use the standard GASM. 1150b1a0a8a3SJed Brown 1151b1a0a8a3SJed Brown Notes: Each processor can have one or more blocks, but a block cannot be shared by more 1152b1a0a8a3SJed Brown than one processor. Defaults to one block per processor. 1153b1a0a8a3SJed Brown 1154b1a0a8a3SJed Brown To set options on the solvers for each block append -sub_ to all the KSP, and PC 1155b1a0a8a3SJed Brown options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly 1156b1a0a8a3SJed Brown 1157f746d493SDmitry Karpeev To set the options on the solvers separate for each block call PCGASMGetSubKSP() 1158b1a0a8a3SJed Brown and set the options directly on the resulting KSP object (you can access its PC 1159b1a0a8a3SJed Brown with KSPGetPC()) 1160b1a0a8a3SJed Brown 1161b1a0a8a3SJed Brown 1162b1a0a8a3SJed Brown Level: beginner 1163b1a0a8a3SJed Brown 1164b1a0a8a3SJed Brown Concepts: additive Schwarz method 1165b1a0a8a3SJed Brown 1166b1a0a8a3SJed Brown References: 1167b1a0a8a3SJed Brown An additive variant of the Schwarz alternating method for the case of many subregions 1168b1a0a8a3SJed Brown M Dryja, OB Widlund - Courant Institute, New York University Technical report 1169b1a0a8a3SJed Brown 1170b1a0a8a3SJed Brown Domain Decompositions: Parallel Multilevel Methods for Elliptic Partial Differential Equations, 1171b1a0a8a3SJed Brown Barry Smith, Petter Bjorstad, and William Gropp, Cambridge University Press, ISBN 0-521-49589-X. 1172b1a0a8a3SJed Brown 1173b1a0a8a3SJed Brown .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 117406b43e7eSDmitry Karpeev PCBJACOBI, PCGASMSetUseTrueLocal(), PCGASMGetSubKSP(), PCGASMSetSubdomains(), 11756a4f0f73SDmitry Karpeev PCGASMSetTotalSubdomains(), PCSetModifySubmatrices(), PCGASMSetOverlap(), PCGASMSetType() 1176b1a0a8a3SJed Brown 1177b1a0a8a3SJed Brown M*/ 1178b1a0a8a3SJed Brown 1179b1a0a8a3SJed Brown EXTERN_C_BEGIN 1180b1a0a8a3SJed Brown #undef __FUNCT__ 1181f746d493SDmitry Karpeev #define __FUNCT__ "PCCreate_GASM" 11827087cfbeSBarry Smith PetscErrorCode PCCreate_GASM(PC pc) 1183b1a0a8a3SJed Brown { 1184b1a0a8a3SJed Brown PetscErrorCode ierr; 1185f746d493SDmitry Karpeev PC_GASM *osm; 1186b1a0a8a3SJed Brown 1187b1a0a8a3SJed Brown PetscFunctionBegin; 1188f746d493SDmitry Karpeev ierr = PetscNewLog(pc,PC_GASM,&osm);CHKERRQ(ierr); 1189ab3e923fSDmitry Karpeev osm->N = PETSC_DECIDE; 119006b43e7eSDmitry Karpeev osm->n = PETSC_DECIDE; 1191ab3e923fSDmitry Karpeev osm->nmax = 0; 1192b1a0a8a3SJed Brown osm->overlap = 1; 1193b1a0a8a3SJed Brown osm->ksp = 0; 11946a4f0f73SDmitry Karpeev osm->gorestriction = 0; 11956a4f0f73SDmitry Karpeev osm->girestriction = 0; 1196ab3e923fSDmitry Karpeev osm->gx = 0; 1197ab3e923fSDmitry Karpeev osm->gy = 0; 1198b1a0a8a3SJed Brown osm->x = 0; 1199b1a0a8a3SJed Brown osm->y = 0; 12006a4f0f73SDmitry Karpeev osm->ois = 0; 12016a4f0f73SDmitry Karpeev osm->iis = 0; 1202b1a0a8a3SJed Brown osm->pmat = 0; 1203f746d493SDmitry Karpeev osm->type = PC_GASM_RESTRICT; 12046a4f0f73SDmitry Karpeev osm->same_subdomain_solvers = PETSC_TRUE; 1205b1a0a8a3SJed Brown osm->sort_indices = PETSC_TRUE; 1206b1a0a8a3SJed Brown 1207b1a0a8a3SJed Brown pc->data = (void*)osm; 1208f746d493SDmitry Karpeev pc->ops->apply = PCApply_GASM; 1209f746d493SDmitry Karpeev pc->ops->applytranspose = PCApplyTranspose_GASM; 1210f746d493SDmitry Karpeev pc->ops->setup = PCSetUp_GASM; 1211a06653b4SBarry Smith pc->ops->reset = PCReset_GASM; 1212f746d493SDmitry Karpeev pc->ops->destroy = PCDestroy_GASM; 1213f746d493SDmitry Karpeev pc->ops->setfromoptions = PCSetFromOptions_GASM; 1214f746d493SDmitry Karpeev pc->ops->setuponblocks = PCSetUpOnBlocks_GASM; 1215f746d493SDmitry Karpeev pc->ops->view = PCView_GASM; 1216b1a0a8a3SJed Brown pc->ops->applyrichardson = 0; 1217b1a0a8a3SJed Brown 121806b43e7eSDmitry Karpeev ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMSetSubdomains_C","PCGASMSetSubdomains_GASM", 121906b43e7eSDmitry Karpeev PCGASMSetSubdomains_GASM);CHKERRQ(ierr); 12206a4f0f73SDmitry Karpeev ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMSetTotalSubdomains_C","PCGASMSetTotalSubdomains_GASM", 12216a4f0f73SDmitry Karpeev PCGASMSetTotalSubdomains_GASM);CHKERRQ(ierr); 1222f746d493SDmitry Karpeev ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMSetOverlap_C","PCGASMSetOverlap_GASM", 1223f746d493SDmitry Karpeev PCGASMSetOverlap_GASM);CHKERRQ(ierr); 1224f746d493SDmitry Karpeev ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMSetType_C","PCGASMSetType_GASM", 1225f746d493SDmitry Karpeev PCGASMSetType_GASM);CHKERRQ(ierr); 1226f746d493SDmitry Karpeev ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMSetSortIndices_C","PCGASMSetSortIndices_GASM", 1227f746d493SDmitry Karpeev PCGASMSetSortIndices_GASM);CHKERRQ(ierr); 1228f746d493SDmitry Karpeev ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMGetSubKSP_C","PCGASMGetSubKSP_GASM", 1229f746d493SDmitry Karpeev PCGASMGetSubKSP_GASM);CHKERRQ(ierr); 1230b1a0a8a3SJed Brown PetscFunctionReturn(0); 1231b1a0a8a3SJed Brown } 1232b1a0a8a3SJed Brown EXTERN_C_END 1233b1a0a8a3SJed Brown 1234b1a0a8a3SJed Brown 1235b1a0a8a3SJed Brown #undef __FUNCT__ 123606b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMCreateLocalSubdomains" 1237b1a0a8a3SJed Brown /*@C 123806b43e7eSDmitry Karpeev PCGASMCreateLocalSubdomains - Creates n local index sets for the overlapping 123906b43e7eSDmitry Karpeev Schwarz preconditioner for a any problem based on its matrix. 1240b1a0a8a3SJed Brown 1241b1a0a8a3SJed Brown Collective 1242b1a0a8a3SJed Brown 1243b1a0a8a3SJed Brown Input Parameters: 1244b1a0a8a3SJed Brown + A - The global matrix operator 12456a4f0f73SDmitry Karpeev . overlap - amount of overlap in outer subdomains 124606b43e7eSDmitry Karpeev - n - the number of local subdomains 1247b1a0a8a3SJed Brown 1248b1a0a8a3SJed Brown Output Parameters: 12496a4f0f73SDmitry Karpeev + iis - the array of index sets defining the local inner subdomains (on which the correction is applied) 12506a4f0f73SDmitry Karpeev - ois - the array of index sets defining the local outer subdomains (on which the residual is computed) 1251b1a0a8a3SJed Brown 1252b1a0a8a3SJed Brown Level: advanced 1253b1a0a8a3SJed Brown 12546a4f0f73SDmitry Karpeev Note: this generates n nonoverlapping local inner subdomains on PETSC_COMM_SELF; 12556a4f0f73SDmitry Karpeev PCGASM will generate the overlap from these if you use them in PCGASMSetSubdomains() and set a 12566a4f0f73SDmitry Karpeev nonzero overlap with PCGASMSetOverlap() 1257b1a0a8a3SJed Brown 1258b1a0a8a3SJed Brown In the Fortran version you must provide the array outis[] already allocated of length n. 1259b1a0a8a3SJed Brown 1260f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, create, subdomains, unstructured grid 1261b1a0a8a3SJed Brown 126206b43e7eSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMDestroySubdomains() 1263b1a0a8a3SJed Brown @*/ 12646a4f0f73SDmitry Karpeev PetscErrorCode PCGASMCreateLocalSubdomains(Mat A, PetscInt overlap, PetscInt n, IS* iis[], IS* ois[]) 1265b1a0a8a3SJed Brown { 1266b1a0a8a3SJed Brown MatPartitioning mpart; 1267b1a0a8a3SJed Brown const char *prefix; 126811bd1e4dSLisandro Dalcin PetscErrorCode (*f)(Mat,MatReuse,Mat*); 1269b1a0a8a3SJed Brown PetscMPIInt size; 1270b1a0a8a3SJed Brown PetscInt i,j,rstart,rend,bs; 127111bd1e4dSLisandro Dalcin PetscBool isbaij = PETSC_FALSE,foundpart = PETSC_FALSE; 1272b1a0a8a3SJed Brown Mat Ad = PETSC_NULL, adj; 1273b1a0a8a3SJed Brown IS ispart,isnumb,*is; 1274b1a0a8a3SJed Brown PetscErrorCode ierr; 1275b1a0a8a3SJed Brown 1276b1a0a8a3SJed Brown PetscFunctionBegin; 1277b1a0a8a3SJed Brown PetscValidHeaderSpecific(A,MAT_CLASSID,1); 12786a4f0f73SDmitry Karpeev PetscValidPointer(iis,4); 1279b1a0a8a3SJed Brown if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of local blocks must be > 0, n = %D",n); 1280b1a0a8a3SJed Brown 1281b1a0a8a3SJed Brown /* Get prefix, row distribution, and block size */ 1282b1a0a8a3SJed Brown ierr = MatGetOptionsPrefix(A,&prefix);CHKERRQ(ierr); 1283b1a0a8a3SJed Brown ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 1284b1a0a8a3SJed Brown ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr); 1285b1a0a8a3SJed 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); 1286b1a0a8a3SJed Brown 1287b1a0a8a3SJed Brown /* Get diagonal block from matrix if possible */ 1288b1a0a8a3SJed Brown ierr = MPI_Comm_size(((PetscObject)A)->comm,&size);CHKERRQ(ierr); 1289b1a0a8a3SJed Brown ierr = PetscObjectQueryFunction((PetscObject)A,"MatGetDiagonalBlock_C",(void (**)(void))&f);CHKERRQ(ierr); 1290b1a0a8a3SJed Brown if (f) { 129111bd1e4dSLisandro Dalcin ierr = MatGetDiagonalBlock(A,&Ad);CHKERRQ(ierr); 1292b1a0a8a3SJed Brown } else if (size == 1) { 129311bd1e4dSLisandro Dalcin Ad = A; 1294b1a0a8a3SJed Brown } 1295b1a0a8a3SJed Brown if (Ad) { 1296251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)Ad,MATSEQBAIJ,&isbaij);CHKERRQ(ierr); 1297251f4c67SDmitry Karpeev if (!isbaij) {ierr = PetscObjectTypeCompare((PetscObject)Ad,MATSEQSBAIJ,&isbaij);CHKERRQ(ierr);} 1298b1a0a8a3SJed Brown } 1299b1a0a8a3SJed Brown if (Ad && n > 1) { 1300b1a0a8a3SJed Brown PetscBool match,done; 1301b1a0a8a3SJed Brown /* Try to setup a good matrix partitioning if available */ 1302b1a0a8a3SJed Brown ierr = MatPartitioningCreate(PETSC_COMM_SELF,&mpart);CHKERRQ(ierr); 1303b1a0a8a3SJed Brown ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr); 1304b1a0a8a3SJed Brown ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr); 1305251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGCURRENT,&match);CHKERRQ(ierr); 1306b1a0a8a3SJed Brown if (!match) { 1307251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGSQUARE,&match);CHKERRQ(ierr); 1308b1a0a8a3SJed Brown } 1309b1a0a8a3SJed Brown if (!match) { /* assume a "good" partitioner is available */ 13101a83f524SJed Brown PetscInt na; 13111a83f524SJed Brown const PetscInt *ia,*ja; 1312b1a0a8a3SJed Brown ierr = MatGetRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr); 1313b1a0a8a3SJed Brown if (done) { 1314b1a0a8a3SJed Brown /* Build adjacency matrix by hand. Unfortunately a call to 1315b1a0a8a3SJed Brown MatConvert(Ad,MATMPIADJ,MAT_INITIAL_MATRIX,&adj) will 1316b1a0a8a3SJed Brown remove the block-aij structure and we cannot expect 1317b1a0a8a3SJed Brown MatPartitioning to split vertices as we need */ 13181a83f524SJed Brown PetscInt i,j,len,nnz,cnt,*iia=0,*jja=0; 13191a83f524SJed Brown const PetscInt *row; 1320b1a0a8a3SJed Brown nnz = 0; 1321b1a0a8a3SJed Brown for (i=0; i<na; i++) { /* count number of nonzeros */ 1322b1a0a8a3SJed Brown len = ia[i+1] - ia[i]; 1323b1a0a8a3SJed Brown row = ja + ia[i]; 1324b1a0a8a3SJed Brown for (j=0; j<len; j++) { 1325b1a0a8a3SJed Brown if (row[j] == i) { /* don't count diagonal */ 1326b1a0a8a3SJed Brown len--; break; 1327b1a0a8a3SJed Brown } 1328b1a0a8a3SJed Brown } 1329b1a0a8a3SJed Brown nnz += len; 1330b1a0a8a3SJed Brown } 1331b1a0a8a3SJed Brown ierr = PetscMalloc((na+1)*sizeof(PetscInt),&iia);CHKERRQ(ierr); 1332b1a0a8a3SJed Brown ierr = PetscMalloc((nnz)*sizeof(PetscInt),&jja);CHKERRQ(ierr); 1333b1a0a8a3SJed Brown nnz = 0; 1334b1a0a8a3SJed Brown iia[0] = 0; 1335b1a0a8a3SJed Brown for (i=0; i<na; i++) { /* fill adjacency */ 1336b1a0a8a3SJed Brown cnt = 0; 1337b1a0a8a3SJed Brown len = ia[i+1] - ia[i]; 1338b1a0a8a3SJed Brown row = ja + ia[i]; 1339b1a0a8a3SJed Brown for (j=0; j<len; j++) { 1340b1a0a8a3SJed Brown if (row[j] != i) { /* if not diagonal */ 1341b1a0a8a3SJed Brown jja[nnz+cnt++] = row[j]; 1342b1a0a8a3SJed Brown } 1343b1a0a8a3SJed Brown } 1344b1a0a8a3SJed Brown nnz += cnt; 1345b1a0a8a3SJed Brown iia[i+1] = nnz; 1346b1a0a8a3SJed Brown } 1347b1a0a8a3SJed Brown /* Partitioning of the adjacency matrix */ 1348b1a0a8a3SJed Brown ierr = MatCreateMPIAdj(PETSC_COMM_SELF,na,na,iia,jja,PETSC_NULL,&adj);CHKERRQ(ierr); 1349b1a0a8a3SJed Brown ierr = MatPartitioningSetAdjacency(mpart,adj);CHKERRQ(ierr); 1350b1a0a8a3SJed Brown ierr = MatPartitioningSetNParts(mpart,n);CHKERRQ(ierr); 1351b1a0a8a3SJed Brown ierr = MatPartitioningApply(mpart,&ispart);CHKERRQ(ierr); 1352b1a0a8a3SJed Brown ierr = ISPartitioningToNumbering(ispart,&isnumb);CHKERRQ(ierr); 13536bf464f9SBarry Smith ierr = MatDestroy(&adj);CHKERRQ(ierr); 1354b1a0a8a3SJed Brown foundpart = PETSC_TRUE; 1355b1a0a8a3SJed Brown } 1356b1a0a8a3SJed Brown ierr = MatRestoreRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr); 1357b1a0a8a3SJed Brown } 13586bf464f9SBarry Smith ierr = MatPartitioningDestroy(&mpart);CHKERRQ(ierr); 1359b1a0a8a3SJed Brown } 1360b1a0a8a3SJed Brown ierr = PetscMalloc(n*sizeof(IS),&is);CHKERRQ(ierr); 1361b1a0a8a3SJed Brown if (!foundpart) { 1362b1a0a8a3SJed Brown 1363b1a0a8a3SJed Brown /* Partitioning by contiguous chunks of rows */ 1364b1a0a8a3SJed Brown 1365b1a0a8a3SJed Brown PetscInt mbs = (rend-rstart)/bs; 1366b1a0a8a3SJed Brown PetscInt start = rstart; 1367b1a0a8a3SJed Brown for (i=0; i<n; i++) { 1368b1a0a8a3SJed Brown PetscInt count = (mbs/n + ((mbs % n) > i)) * bs; 1369b1a0a8a3SJed Brown ierr = ISCreateStride(PETSC_COMM_SELF,count,start,1,&is[i]);CHKERRQ(ierr); 1370b1a0a8a3SJed Brown start += count; 1371b1a0a8a3SJed Brown } 1372b1a0a8a3SJed Brown 1373b1a0a8a3SJed Brown } else { 1374b1a0a8a3SJed Brown 1375b1a0a8a3SJed Brown /* Partitioning by adjacency of diagonal block */ 1376b1a0a8a3SJed Brown 1377b1a0a8a3SJed Brown const PetscInt *numbering; 1378b1a0a8a3SJed Brown PetscInt *count,nidx,*indices,*newidx,start=0; 1379b1a0a8a3SJed Brown /* Get node count in each partition */ 1380b1a0a8a3SJed Brown ierr = PetscMalloc(n*sizeof(PetscInt),&count);CHKERRQ(ierr); 1381b1a0a8a3SJed Brown ierr = ISPartitioningCount(ispart,n,count);CHKERRQ(ierr); 1382b1a0a8a3SJed Brown if (isbaij && bs > 1) { /* adjust for the block-aij case */ 1383b1a0a8a3SJed Brown for (i=0; i<n; i++) count[i] *= bs; 1384b1a0a8a3SJed Brown } 1385b1a0a8a3SJed Brown /* Build indices from node numbering */ 1386b1a0a8a3SJed Brown ierr = ISGetLocalSize(isnumb,&nidx);CHKERRQ(ierr); 1387b1a0a8a3SJed Brown ierr = PetscMalloc(nidx*sizeof(PetscInt),&indices);CHKERRQ(ierr); 1388b1a0a8a3SJed Brown for (i=0; i<nidx; i++) indices[i] = i; /* needs to be initialized */ 1389b1a0a8a3SJed Brown ierr = ISGetIndices(isnumb,&numbering);CHKERRQ(ierr); 1390b1a0a8a3SJed Brown ierr = PetscSortIntWithPermutation(nidx,numbering,indices);CHKERRQ(ierr); 1391b1a0a8a3SJed Brown ierr = ISRestoreIndices(isnumb,&numbering);CHKERRQ(ierr); 1392b1a0a8a3SJed Brown if (isbaij && bs > 1) { /* adjust for the block-aij case */ 1393b1a0a8a3SJed Brown ierr = PetscMalloc(nidx*bs*sizeof(PetscInt),&newidx);CHKERRQ(ierr); 1394b1a0a8a3SJed Brown for (i=0; i<nidx; i++) 1395b1a0a8a3SJed Brown for (j=0; j<bs; j++) 1396b1a0a8a3SJed Brown newidx[i*bs+j] = indices[i]*bs + j; 1397b1a0a8a3SJed Brown ierr = PetscFree(indices);CHKERRQ(ierr); 1398b1a0a8a3SJed Brown nidx *= bs; 1399b1a0a8a3SJed Brown indices = newidx; 1400b1a0a8a3SJed Brown } 1401b1a0a8a3SJed Brown /* Shift to get global indices */ 1402b1a0a8a3SJed Brown for (i=0; i<nidx; i++) indices[i] += rstart; 1403b1a0a8a3SJed Brown 1404b1a0a8a3SJed Brown /* Build the index sets for each block */ 1405b1a0a8a3SJed Brown for (i=0; i<n; i++) { 1406b1a0a8a3SJed Brown ierr = ISCreateGeneral(PETSC_COMM_SELF,count[i],&indices[start],PETSC_COPY_VALUES,&is[i]);CHKERRQ(ierr); 1407b1a0a8a3SJed Brown ierr = ISSort(is[i]);CHKERRQ(ierr); 1408b1a0a8a3SJed Brown start += count[i]; 1409b1a0a8a3SJed Brown } 1410b1a0a8a3SJed Brown 1411b1a0a8a3SJed Brown ierr = PetscFree(count); 1412b1a0a8a3SJed Brown ierr = PetscFree(indices); 1413fcfd50ebSBarry Smith ierr = ISDestroy(&isnumb);CHKERRQ(ierr); 1414fcfd50ebSBarry Smith ierr = ISDestroy(&ispart);CHKERRQ(ierr); 1415b1a0a8a3SJed Brown } 14166a4f0f73SDmitry Karpeev *iis = is; 14176a4f0f73SDmitry Karpeev if (!ois) PetscFunctionReturn(0); 14186a4f0f73SDmitry Karpeev /* 14196a4f0f73SDmitry Karpeev Initially make outer subdomains the same as inner subdomains. If nonzero additional overlap 14206a4f0f73SDmitry Karpeev has been requested, copy the inner subdomains over so they can be modified. 14216a4f0f73SDmitry Karpeev */ 14226a4f0f73SDmitry Karpeev ierr = PetscMalloc(n*sizeof(IS),ois);CHKERRQ(ierr); 14236a4f0f73SDmitry Karpeev for (i=0; i<n; ++i) { 14246a4f0f73SDmitry Karpeev if (overlap > 0) { /* With positive overlap, (*iis)[i] will be modified */ 14256a4f0f73SDmitry Karpeev ierr = ISDuplicate((*iis)[i],(*ois)+i);CHKERRQ(ierr); 14266a4f0f73SDmitry Karpeev ierr = ISCopy((*iis)[i],(*ois)[i]);CHKERRQ(ierr); 14276a4f0f73SDmitry Karpeev } else { 14286a4f0f73SDmitry Karpeev ierr = PetscObjectReference((PetscObject)(*iis)[i]);CHKERRQ(ierr); 14296a4f0f73SDmitry Karpeev (*ois)[i] = (*iis)[i]; 14306a4f0f73SDmitry Karpeev } 14316a4f0f73SDmitry Karpeev } 14326a4f0f73SDmitry Karpeev if (overlap > 0) { 14336a4f0f73SDmitry Karpeev /* Extend the "overlapping" regions by a number of steps */ 14346a4f0f73SDmitry Karpeev ierr = MatIncreaseOverlap(A,n,*ois,overlap);CHKERRQ(ierr); 14356a4f0f73SDmitry Karpeev } 1436b1a0a8a3SJed Brown PetscFunctionReturn(0); 1437b1a0a8a3SJed Brown } 1438b1a0a8a3SJed Brown 1439b1a0a8a3SJed Brown #undef __FUNCT__ 1440f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMDestroySubdomains" 1441b1a0a8a3SJed Brown /*@C 1442f746d493SDmitry Karpeev PCGASMDestroySubdomains - Destroys the index sets created with 144306b43e7eSDmitry Karpeev PCGASMCreateLocalSubdomains() or PCGASMCreateSubdomains2D. Should be 144406b43e7eSDmitry Karpeev called after setting subdomains with PCGASMSetSubdomains(). 1445b1a0a8a3SJed Brown 1446b1a0a8a3SJed Brown Collective 1447b1a0a8a3SJed Brown 1448b1a0a8a3SJed Brown Input Parameters: 1449b1a0a8a3SJed Brown + n - the number of index sets 14506a4f0f73SDmitry Karpeev . iis - the array of inner subdomains, 14516a4f0f73SDmitry Karpeev - ois - the array of outer subdomains, can be PETSC_NULL 1452b1a0a8a3SJed Brown 14536a4f0f73SDmitry Karpeev Level: intermediate 14546a4f0f73SDmitry Karpeev 14556a4f0f73SDmitry Karpeev Notes: this is merely a convenience subroutine that walks each list, 14566a4f0f73SDmitry Karpeev destroys each IS on the list, and then frees the list. 1457b1a0a8a3SJed Brown 1458f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, create, subdomains, unstructured grid 1459b1a0a8a3SJed Brown 146006b43e7eSDmitry Karpeev .seealso: PCGASMCreateLocalSubdomains(), PCGASMSetSubdomains() 1461b1a0a8a3SJed Brown @*/ 14626a4f0f73SDmitry Karpeev PetscErrorCode PCGASMDestroySubdomains(PetscInt n, IS iis[], IS ois[]) 1463b1a0a8a3SJed Brown { 1464b1a0a8a3SJed Brown PetscInt i; 1465b1a0a8a3SJed Brown PetscErrorCode ierr; 1466b1a0a8a3SJed Brown PetscFunctionBegin; 1467a35b7b57SDmitry Karpeev if (n <= 0) PetscFunctionReturn(0); 1468a35b7b57SDmitry Karpeev if (iis) { 14696a4f0f73SDmitry Karpeev PetscValidPointer(iis,2); 1470a35b7b57SDmitry Karpeev for (i=0; i<n; i++) { 1471a35b7b57SDmitry Karpeev ierr = ISDestroy(&iis[i]);CHKERRQ(ierr); 1472a35b7b57SDmitry Karpeev } 14736a4f0f73SDmitry Karpeev ierr = PetscFree(iis);CHKERRQ(ierr); 1474a35b7b57SDmitry Karpeev } 14756a4f0f73SDmitry Karpeev if (ois) { 1476a35b7b57SDmitry Karpeev for (i=0; i<n; i++) { 1477a35b7b57SDmitry Karpeev ierr = ISDestroy(&ois[i]);CHKERRQ(ierr); 1478a35b7b57SDmitry Karpeev } 14796a4f0f73SDmitry Karpeev ierr = PetscFree(ois);CHKERRQ(ierr); 1480b1a0a8a3SJed Brown } 1481b1a0a8a3SJed Brown PetscFunctionReturn(0); 1482b1a0a8a3SJed Brown } 1483b1a0a8a3SJed Brown 1484af538404SDmitry Karpeev 1485af538404SDmitry Karpeev #define PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,xleft_loc,ylow_loc,xright_loc,yhigh_loc,n) \ 1486af538404SDmitry Karpeev { \ 1487af538404SDmitry Karpeev PetscInt first_row = first/M, last_row = last/M+1; \ 1488af538404SDmitry Karpeev /* \ 1489af538404SDmitry Karpeev Compute ylow_loc and yhigh_loc so that (ylow_loc,xleft) and (yhigh_loc,xright) are the corners \ 1490af538404SDmitry Karpeev of the bounding box of the intersection of the subdomain with the local ownership range (local \ 1491af538404SDmitry Karpeev subdomain). \ 1492af538404SDmitry Karpeev Also compute xleft_loc and xright_loc as the lower and upper bounds on the first and last rows \ 1493af538404SDmitry Karpeev of the intersection. \ 1494af538404SDmitry Karpeev */ \ 1495af538404SDmitry Karpeev /* ylow_loc is the grid row containing the first element of the local sumbdomain */ \ 1496eec7e3faSDmitry Karpeev *ylow_loc = PetscMax(first_row,ylow); \ 1497af538404SDmitry Karpeev /* xleft_loc is the offset of first element of the local subdomain within its grid row (might actually be outside the local subdomain) */ \ 1498eec7e3faSDmitry Karpeev *xleft_loc = *ylow_loc==first_row?PetscMax(first%M,xleft):xleft; \ 1499af538404SDmitry Karpeev /* yhigh_loc is the grid row above the last local subdomain element */ \ 1500eec7e3faSDmitry Karpeev *yhigh_loc = PetscMin(last_row,yhigh); \ 1501af538404SDmitry Karpeev /* xright is the offset of the end of the local subdomain within its grid row (might actually be outside the local subdomain) */ \ 1502eec7e3faSDmitry Karpeev *xright_loc = *yhigh_loc==last_row?PetscMin(xright,last%M):xright; \ 1503af538404SDmitry Karpeev /* Now compute the size of the local subdomain n. */ \ 1504c3518bceSDmitry Karpeev *n = 0; \ 1505eec7e3faSDmitry Karpeev if (*ylow_loc < *yhigh_loc) { \ 1506af538404SDmitry Karpeev PetscInt width = xright-xleft; \ 1507eec7e3faSDmitry Karpeev *n += width*(*yhigh_loc-*ylow_loc-1); \ 1508eec7e3faSDmitry Karpeev *n += PetscMin(PetscMax(*xright_loc-xleft,0),width); \ 1509eec7e3faSDmitry Karpeev *n -= PetscMin(PetscMax(*xleft_loc-xleft,0), width); \ 1510af538404SDmitry Karpeev }\ 1511af538404SDmitry Karpeev } 1512af538404SDmitry Karpeev 1513c3518bceSDmitry Karpeev 1514eec7e3faSDmitry Karpeev 1515eec7e3faSDmitry Karpeev #undef __FUNCT__ 1516f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMCreateSubdomains2D" 1517b1a0a8a3SJed Brown /*@ 1518f746d493SDmitry Karpeev PCGASMCreateSubdomains2D - Creates the index sets for the overlapping Schwarz 1519b1a0a8a3SJed Brown preconditioner for a two-dimensional problem on a regular grid. 1520b1a0a8a3SJed Brown 1521af538404SDmitry Karpeev Collective 1522b1a0a8a3SJed Brown 1523b1a0a8a3SJed Brown Input Parameters: 152406b43e7eSDmitry Karpeev + M, N - the global number of grid points in the x and y directions 1525af538404SDmitry Karpeev . Mdomains, Ndomains - the global number of subdomains in the x and y directions 1526b1a0a8a3SJed Brown . dof - degrees of freedom per node 1527b1a0a8a3SJed Brown - overlap - overlap in mesh lines 1528b1a0a8a3SJed Brown 1529b1a0a8a3SJed Brown Output Parameters: 1530af538404SDmitry Karpeev + Nsub - the number of local subdomains created 15316a4f0f73SDmitry Karpeev . iis - array of index sets defining inner (nonoverlapping) subdomains 15326a4f0f73SDmitry Karpeev - ois - array of index sets defining outer (overlapping, if overlap > 0) subdomains 1533b1a0a8a3SJed Brown 1534b1a0a8a3SJed Brown 1535b1a0a8a3SJed Brown Level: advanced 1536b1a0a8a3SJed Brown 1537f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, create, subdomains, 2D, regular grid 1538b1a0a8a3SJed Brown 15396a4f0f73SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetSubdomains(), PCGASMGetSubKSP(), 1540f746d493SDmitry Karpeev PCGASMSetOverlap() 1541b1a0a8a3SJed Brown @*/ 15426a4f0f73SDmitry Karpeev PetscErrorCode PCGASMCreateSubdomains2D(PC pc, PetscInt M,PetscInt N,PetscInt Mdomains,PetscInt Ndomains,PetscInt dof,PetscInt overlap, PetscInt *nsub,IS **iis,IS **ois) 1543b1a0a8a3SJed Brown { 1544b1a0a8a3SJed Brown PetscErrorCode ierr; 15452bbb417fSDmitry Karpeev PetscMPIInt size, rank; 15462bbb417fSDmitry Karpeev PetscInt i, j; 15472bbb417fSDmitry Karpeev PetscInt maxheight, maxwidth; 15482bbb417fSDmitry Karpeev PetscInt xstart, xleft, xright, xleft_loc, xright_loc; 15492bbb417fSDmitry Karpeev PetscInt ystart, ylow, yhigh, ylow_loc, yhigh_loc; 1550eec7e3faSDmitry Karpeev PetscInt x[2][2], y[2][2], n[2]; 15512bbb417fSDmitry Karpeev PetscInt first, last; 15522bbb417fSDmitry Karpeev PetscInt nidx, *idx; 15532bbb417fSDmitry Karpeev PetscInt ii,jj,s,q,d; 1554af538404SDmitry Karpeev PetscInt k,kk; 15552bbb417fSDmitry Karpeev PetscMPIInt color; 15562bbb417fSDmitry Karpeev MPI_Comm comm, subcomm; 15576a4f0f73SDmitry Karpeev IS **xis = 0, **is = ois, **is_local = iis; 1558d34fcf5fSBarry Smith 1559b1a0a8a3SJed Brown PetscFunctionBegin; 15602bbb417fSDmitry Karpeev ierr = PetscObjectGetComm((PetscObject)pc, &comm);CHKERRQ(ierr); 15612bbb417fSDmitry Karpeev ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 15622bbb417fSDmitry Karpeev ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 15632bbb417fSDmitry Karpeev ierr = MatGetOwnershipRange(pc->pmat, &first, &last);CHKERRQ(ierr); 1564d34fcf5fSBarry 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) " 15652bbb417fSDmitry Karpeev "does not respect the number of degrees of freedom per grid point %D", first, last, dof); 1566d34fcf5fSBarry Smith 1567af538404SDmitry Karpeev /* Determine the number of domains with nonzero intersections with the local ownership range. */ 15682bbb417fSDmitry Karpeev s = 0; 1569af538404SDmitry Karpeev ystart = 0; 1570af538404SDmitry Karpeev for (j=0; j<Ndomains; ++j) { 1571af538404SDmitry Karpeev maxheight = N/Ndomains + ((N % Ndomains) > j); /* Maximal height of subdomain */ 1572af538404SDmitry 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); 1573eec7e3faSDmitry Karpeev /* Vertical domain limits with an overlap. */ 1574eec7e3faSDmitry Karpeev ylow = PetscMax(ystart - overlap,0); 1575eec7e3faSDmitry Karpeev yhigh = PetscMin(ystart + maxheight + overlap,N); 1576b1a0a8a3SJed Brown xstart = 0; 1577af538404SDmitry Karpeev for (i=0; i<Mdomains; ++i) { 1578af538404SDmitry Karpeev maxwidth = M/Mdomains + ((M % Mdomains) > i); /* Maximal width of subdomain */ 1579af538404SDmitry 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); 1580eec7e3faSDmitry Karpeev /* Horizontal domain limits with an overlap. */ 1581eec7e3faSDmitry Karpeev xleft = PetscMax(xstart - overlap,0); 1582eec7e3faSDmitry Karpeev xright = PetscMin(xstart + maxwidth + overlap,M); 1583af538404SDmitry Karpeev /* 1584af538404SDmitry Karpeev Determine whether this subdomain intersects this processor's ownership range of pc->pmat. 1585af538404SDmitry Karpeev */ 1586c3518bceSDmitry Karpeev PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,(&xleft_loc),(&ylow_loc),(&xright_loc),(&yhigh_loc),(&nidx)); 1587af538404SDmitry Karpeev if (nidx) { 1588af538404SDmitry Karpeev ++s; 1589af538404SDmitry Karpeev } 1590af538404SDmitry Karpeev xstart += maxwidth; 1591af538404SDmitry Karpeev }/* for (i = 0; i < Mdomains; ++i) */ 1592af538404SDmitry Karpeev ystart += maxheight; 1593af538404SDmitry Karpeev }/* for (j = 0; j < Ndomains; ++j) */ 1594af538404SDmitry Karpeev /* Now we can allocate the necessary number of ISs. */ 1595af538404SDmitry Karpeev *nsub = s; 1596af538404SDmitry Karpeev ierr = PetscMalloc((*nsub)*sizeof(IS*),is);CHKERRQ(ierr); 1597af538404SDmitry Karpeev ierr = PetscMalloc((*nsub)*sizeof(IS*),is_local);CHKERRQ(ierr); 1598af538404SDmitry Karpeev s = 0; 1599af538404SDmitry Karpeev ystart = 0; 1600af538404SDmitry Karpeev for (j=0; j<Ndomains; ++j) { 1601af538404SDmitry Karpeev maxheight = N/Ndomains + ((N % Ndomains) > j); /* Maximal height of subdomain */ 1602af538404SDmitry 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); 1603af538404SDmitry Karpeev /* Vertical domain limits with an overlap. */ 1604eec7e3faSDmitry Karpeev y[0][0] = PetscMax(ystart - overlap,0); 1605eec7e3faSDmitry Karpeev y[0][1] = PetscMin(ystart + maxheight + overlap,N); 1606eec7e3faSDmitry Karpeev /* Vertical domain limits without an overlap. */ 1607eec7e3faSDmitry Karpeev y[1][0] = ystart; 1608eec7e3faSDmitry Karpeev y[1][1] = PetscMin(ystart + maxheight,N); 1609eec7e3faSDmitry Karpeev xstart = 0; 1610af538404SDmitry Karpeev for (i=0; i<Mdomains; ++i) { 1611af538404SDmitry Karpeev maxwidth = M/Mdomains + ((M % Mdomains) > i); /* Maximal width of subdomain */ 1612af538404SDmitry 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); 1613af538404SDmitry Karpeev /* Horizontal domain limits with an overlap. */ 1614eec7e3faSDmitry Karpeev x[0][0] = PetscMax(xstart - overlap,0); 1615eec7e3faSDmitry Karpeev x[0][1] = PetscMin(xstart + maxwidth + overlap,M); 1616eec7e3faSDmitry Karpeev /* Horizontal domain limits without an overlap. */ 1617eec7e3faSDmitry Karpeev x[1][0] = xstart; 1618eec7e3faSDmitry Karpeev x[1][1] = PetscMin(xstart+maxwidth,M); 16192bbb417fSDmitry Karpeev /* 16202bbb417fSDmitry Karpeev Determine whether this domain intersects this processor's ownership range of pc->pmat. 16212bbb417fSDmitry Karpeev Do this twice: first for the domains with overlaps, and once without. 16222bbb417fSDmitry Karpeev During the first pass create the subcommunicators, and use them on the second pass as well. 16232bbb417fSDmitry Karpeev */ 16242bbb417fSDmitry Karpeev for (q = 0; q < 2; ++q) { 16252bbb417fSDmitry Karpeev /* 16262bbb417fSDmitry Karpeev domain limits, (xleft, xright) and (ylow, yheigh) are adjusted 16272bbb417fSDmitry Karpeev according to whether the domain with an overlap or without is considered. 16282bbb417fSDmitry Karpeev */ 16292bbb417fSDmitry Karpeev xleft = x[q][0]; xright = x[q][1]; 16302bbb417fSDmitry Karpeev ylow = y[q][0]; yhigh = y[q][1]; 1631c3518bceSDmitry Karpeev PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,(&xleft_loc),(&ylow_loc),(&xright_loc),(&yhigh_loc),(&nidx)); 1632af538404SDmitry Karpeev nidx *= dof; 1633eec7e3faSDmitry Karpeev n[q] = nidx; 16342bbb417fSDmitry Karpeev /* 1635eec7e3faSDmitry Karpeev Based on the counted number of indices in the local domain *with an overlap*, 1636af538404SDmitry Karpeev construct a subcommunicator of all the processors supporting this domain. 1637eec7e3faSDmitry Karpeev Observe that a domain with an overlap might have nontrivial local support, 1638eec7e3faSDmitry Karpeev while the domain without an overlap might not. Hence, the decision to participate 1639eec7e3faSDmitry Karpeev in the subcommunicator must be based on the domain with an overlap. 16402bbb417fSDmitry Karpeev */ 16412bbb417fSDmitry Karpeev if (q == 0) { 16422bbb417fSDmitry Karpeev if (nidx) { 16432bbb417fSDmitry Karpeev color = 1; 1644d34fcf5fSBarry Smith } else { 16452bbb417fSDmitry Karpeev color = MPI_UNDEFINED; 16462bbb417fSDmitry Karpeev } 16472bbb417fSDmitry Karpeev ierr = MPI_Comm_split(comm, color, rank, &subcomm);CHKERRQ(ierr); 16482bbb417fSDmitry Karpeev } 1649af538404SDmitry Karpeev /* 1650eec7e3faSDmitry Karpeev Proceed only if the number of local indices *with an overlap* is nonzero. 1651af538404SDmitry Karpeev */ 1652eec7e3faSDmitry Karpeev if (n[0]) { 1653af538404SDmitry Karpeev if (q == 0) { 16546a4f0f73SDmitry Karpeev xis = is; 1655af538404SDmitry Karpeev } 1656af538404SDmitry Karpeev if (q == 1) { 1657af538404SDmitry Karpeev /* 1658eec7e3faSDmitry Karpeev The IS for the no-overlap subdomain shares a communicator with the overlapping domain. 1659af538404SDmitry Karpeev Moreover, if the overlap is zero, the two ISs are identical. 1660af538404SDmitry Karpeev */ 1661b1a0a8a3SJed Brown if (overlap == 0) { 1662eec7e3faSDmitry Karpeev (*is_local)[s] = (*is)[s]; 16632bbb417fSDmitry Karpeev ierr = PetscObjectReference((PetscObject)(*is)[s]);CHKERRQ(ierr); 16642bbb417fSDmitry Karpeev continue; 1665d34fcf5fSBarry Smith } else { 16666a4f0f73SDmitry Karpeev xis = is_local; 1667eec7e3faSDmitry Karpeev subcomm = ((PetscObject)(*is)[s])->comm; 16682bbb417fSDmitry Karpeev } 1669af538404SDmitry Karpeev }/* if (q == 1) */ 1670eec7e3faSDmitry Karpeev idx = PETSC_NULL; 16712bbb417fSDmitry Karpeev ierr = PetscMalloc(nidx*sizeof(PetscInt),&idx);CHKERRQ(ierr); 1672eec7e3faSDmitry Karpeev if (nidx) { 16732bbb417fSDmitry Karpeev k = 0; 16742bbb417fSDmitry Karpeev for (jj=ylow_loc; jj<yhigh_loc; ++jj) { 1675af538404SDmitry Karpeev PetscInt x0 = (jj==ylow_loc)?xleft_loc:xleft; 1676af538404SDmitry Karpeev PetscInt x1 = (jj==yhigh_loc-1)?xright_loc:xright; 1677af538404SDmitry Karpeev kk = dof*(M*jj + x0); 16782bbb417fSDmitry Karpeev for (ii=x0; ii<x1; ++ii) { 16792bbb417fSDmitry Karpeev for (d = 0; d < dof; ++d) { 16802bbb417fSDmitry Karpeev idx[k++] = kk++; 1681b1a0a8a3SJed Brown } 1682b1a0a8a3SJed Brown } 1683b1a0a8a3SJed Brown } 1684eec7e3faSDmitry Karpeev } 16856a4f0f73SDmitry Karpeev ierr = ISCreateGeneral(subcomm,nidx,idx,PETSC_OWN_POINTER,(*xis)+s);CHKERRQ(ierr); 1686eec7e3faSDmitry Karpeev }/* if (n[0]) */ 16872bbb417fSDmitry Karpeev }/* for (q = 0; q < 2; ++q) */ 1688eec7e3faSDmitry Karpeev if (n[0]) { 16892bbb417fSDmitry Karpeev ++s; 1690b1a0a8a3SJed Brown } 1691af538404SDmitry Karpeev xstart += maxwidth; 1692af538404SDmitry Karpeev }/* for (i = 0; i < Mdomains; ++i) */ 16932bbb417fSDmitry Karpeev ystart += maxheight; 1694af538404SDmitry Karpeev }/* for (j = 0; j < Ndomains; ++j) */ 1695b1a0a8a3SJed Brown PetscFunctionReturn(0); 1696b1a0a8a3SJed Brown } 1697b1a0a8a3SJed Brown 1698b1a0a8a3SJed Brown #undef __FUNCT__ 169906b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMGetSubdomains" 1700b1a0a8a3SJed Brown /*@C 170106b43e7eSDmitry Karpeev PCGASMGetSubdomains - Gets the subdomains supported on this processor 170206b43e7eSDmitry Karpeev for the additive Schwarz preconditioner. 1703b1a0a8a3SJed Brown 1704b1a0a8a3SJed Brown Not Collective 1705b1a0a8a3SJed Brown 1706b1a0a8a3SJed Brown Input Parameter: 1707b1a0a8a3SJed Brown . pc - the preconditioner context 1708b1a0a8a3SJed Brown 1709b1a0a8a3SJed Brown Output Parameters: 1710b1a0a8a3SJed Brown + n - the number of subdomains for this processor (default value = 1) 17116a4f0f73SDmitry Karpeev . iis - the index sets that define the inner subdomains (without overlap) supported on this processor (can be PETSC_NULL) 17126a4f0f73SDmitry Karpeev - ois - the index sets that define the outer subdomains (with overlap) supported on this processor (can be PETSC_NULL) 1713b1a0a8a3SJed Brown 1714b1a0a8a3SJed Brown 1715b1a0a8a3SJed Brown Notes: 17166a4f0f73SDmitry Karpeev The user is responsible for destroying the ISs and freeing the returned arrays. 1717b1a0a8a3SJed Brown The IS numbering is in the parallel, global numbering of the vector. 1718b1a0a8a3SJed Brown 1719b1a0a8a3SJed Brown Level: advanced 1720b1a0a8a3SJed Brown 172106b43e7eSDmitry Karpeev .keywords: PC, GASM, get, subdomains, additive Schwarz 1722b1a0a8a3SJed Brown 17236a4f0f73SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetOverlap(), PCGASMGetSubKSP(), 172406b43e7eSDmitry Karpeev PCGASMCreateSubdomains2D(), PCGASMSetSubdomains(), PCGASMGetSubmatrices() 1725b1a0a8a3SJed Brown @*/ 17266a4f0f73SDmitry Karpeev PetscErrorCode PCGASMGetSubdomains(PC pc,PetscInt *n,IS *iis[],IS *ois[]) 1727b1a0a8a3SJed Brown { 1728f746d493SDmitry Karpeev PC_GASM *osm; 1729b1a0a8a3SJed Brown PetscErrorCode ierr; 1730b1a0a8a3SJed Brown PetscBool match; 17316a4f0f73SDmitry Karpeev PetscInt i; 1732b1a0a8a3SJed Brown PetscFunctionBegin; 1733b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1734251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match);CHKERRQ(ierr); 17356a4f0f73SDmitry Karpeev if (!match) 17366a4f0f73SDmitry Karpeev SETERRQ2(((PetscObject)pc)->comm, PETSC_ERR_ARG_WRONG, "Incorrect object type: expected %s, got %s instead", PCGASM, ((PetscObject)pc)->type_name); 1737f746d493SDmitry Karpeev osm = (PC_GASM*)pc->data; 1738ab3e923fSDmitry Karpeev if (n) *n = osm->n; 17396a4f0f73SDmitry Karpeev if (iis) { 17406a4f0f73SDmitry Karpeev ierr = PetscMalloc(osm->n*sizeof(IS), iis);CHKERRQ(ierr); 17416a4f0f73SDmitry Karpeev } 17426a4f0f73SDmitry Karpeev if (ois) { 17436a4f0f73SDmitry Karpeev ierr = PetscMalloc(osm->n*sizeof(IS), ois);CHKERRQ(ierr); 17446a4f0f73SDmitry Karpeev } 17456a4f0f73SDmitry Karpeev if (iis || ois) { 17466a4f0f73SDmitry Karpeev for (i = 0; i < osm->n; ++i) { 17476a4f0f73SDmitry Karpeev if (iis) (*iis)[i] = osm->iis[i]; 17486a4f0f73SDmitry Karpeev if (ois) (*ois)[i] = osm->ois[i]; 17496a4f0f73SDmitry Karpeev } 1750b1a0a8a3SJed Brown } 1751b1a0a8a3SJed Brown PetscFunctionReturn(0); 1752b1a0a8a3SJed Brown } 1753b1a0a8a3SJed Brown 1754b1a0a8a3SJed Brown #undef __FUNCT__ 175506b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMGetSubmatrices" 1756b1a0a8a3SJed Brown /*@C 175706b43e7eSDmitry Karpeev PCGASMGetSubmatrices - Gets the local submatrices (for this processor 1758b1a0a8a3SJed Brown only) for the additive Schwarz preconditioner. 1759b1a0a8a3SJed Brown 1760b1a0a8a3SJed Brown Not Collective 1761b1a0a8a3SJed Brown 1762b1a0a8a3SJed Brown Input Parameter: 1763b1a0a8a3SJed Brown . pc - the preconditioner context 1764b1a0a8a3SJed Brown 1765b1a0a8a3SJed Brown Output Parameters: 1766b1a0a8a3SJed Brown + n - the number of matrices for this processor (default value = 1) 1767b1a0a8a3SJed Brown - mat - the matrices 1768b1a0a8a3SJed Brown 176906b43e7eSDmitry Karpeev Notes: matrices returned by this routine have the same communicators as the index sets (IS) 177006b43e7eSDmitry Karpeev used to define subdomains in PCGASMSetSubdomains(), or PETSC_COMM_SELF, if the 17716a4f0f73SDmitry Karpeev subdomains were defined using PCGASMSetTotalSubdomains(). 1772b1a0a8a3SJed Brown Level: advanced 1773b1a0a8a3SJed Brown 1774f746d493SDmitry Karpeev .keywords: PC, GASM, set, local, subdomains, additive Schwarz, block Jacobi 1775b1a0a8a3SJed Brown 17766a4f0f73SDmitry Karpeev .seealso: PCGASMSetTotalSubdomain(), PCGASMSetOverlap(), PCGASMGetSubKSP(), 177706b43e7eSDmitry Karpeev PCGASMCreateSubdomains2D(), PCGASMSetSubdomains(), PCGASMGetSubdomains() 1778b1a0a8a3SJed Brown @*/ 177906b43e7eSDmitry Karpeev PetscErrorCode PCGASMGetSubmatrices(PC pc,PetscInt *n,Mat *mat[]) 1780b1a0a8a3SJed Brown { 1781f746d493SDmitry Karpeev PC_GASM *osm; 1782b1a0a8a3SJed Brown PetscErrorCode ierr; 1783b1a0a8a3SJed Brown PetscBool match; 1784b1a0a8a3SJed Brown 1785b1a0a8a3SJed Brown PetscFunctionBegin; 1786b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1787b1a0a8a3SJed Brown PetscValidIntPointer(n,2); 1788b1a0a8a3SJed Brown if (mat) PetscValidPointer(mat,3); 1789b1a0a8a3SJed Brown if (!pc->setupcalled) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must call after KSPSetUP() or PCSetUp()."); 1790251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match);CHKERRQ(ierr); 179106b43e7eSDmitry Karpeev if (!match) SETERRQ2(((PetscObject)pc)->comm, PETSC_ERR_ARG_WRONG, "Expected %s, got %s instead", PCGASM, ((PetscObject)pc)->type_name); 1792f746d493SDmitry Karpeev osm = (PC_GASM*)pc->data; 1793ab3e923fSDmitry Karpeev if (n) *n = osm->n; 1794b1a0a8a3SJed Brown if (mat) *mat = osm->pmat; 179506b43e7eSDmitry Karpeev 1796b1a0a8a3SJed Brown PetscFunctionReturn(0); 1797b1a0a8a3SJed Brown } 1798