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); 139*7105d80bSDmitry 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. */ 275fdc48646SDmitry Karpeev DM *domain_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; 298fdc48646SDmitry Karpeev PetscInt num_domains, d; 299fdc48646SDmitry Karpeev char **domain_names; 300fdc48646SDmitry Karpeev IS *domain_is; 301fdc48646SDmitry Karpeev /* Allow the user to request a decomposition DM by name */ 302fdc48646SDmitry Karpeev ierr = PetscStrncpy(ddm_name, "", 1024); CHKERRQ(ierr); 30306b43e7eSDmitry Karpeev ierr = PetscOptionsString("-pc_gasm_decomposition","Name of the DM defining the composition", "PCSetDM", ddm_name, ddm_name,1024,&flg); CHKERRQ(ierr); 304fdc48646SDmitry Karpeev if(flg) { 305fdc48646SDmitry Karpeev ierr = DMCreateDecompositionDM(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 } 312fdc48646SDmitry Karpeev ierr = DMCreateDecomposition(pc->dm, &num_domains, &domain_names, &domain_is, &domain_dm); CHKERRQ(ierr); 31369ca1f37SDmitry Karpeev if(num_domains) { 31406b43e7eSDmitry Karpeev ierr = PCGASMSetSubdomains(pc, num_domains, domain_is, PETSC_NULL);CHKERRQ(ierr); 31569ca1f37SDmitry Karpeev } 316fdc48646SDmitry Karpeev for(d = 0; d < num_domains; ++d) { 317fdc48646SDmitry Karpeev ierr = PetscFree(domain_names[d]); CHKERRQ(ierr); 318fdc48646SDmitry Karpeev ierr = ISDestroy(&domain_is[d]); CHKERRQ(ierr); 319fdc48646SDmitry Karpeev } 320fdc48646SDmitry Karpeev ierr = PetscFree(domain_names);CHKERRQ(ierr); 321fdc48646SDmitry Karpeev ierr = PetscFree(domain_is);CHKERRQ(ierr); 322fdc48646SDmitry Karpeev } 32306b43e7eSDmitry Karpeev if (osm->n == PETSC_DECIDE) { /* still no subdomains; use one per processor */ 32444fe18e5SDmitry Karpeev osm->nmax = osm->n = 1; 325b1a0a8a3SJed Brown ierr = MPI_Comm_size(((PetscObject)pc)->comm,&size);CHKERRQ(ierr); 326f746d493SDmitry Karpeev osm->N = size; 327fdc48646SDmitry Karpeev } 32806b43e7eSDmitry Karpeev } 32906b43e7eSDmitry Karpeev if (osm->N == PETSC_DECIDE) { 330b1a0a8a3SJed Brown PetscInt inwork[2], outwork[2]; 331f746d493SDmitry Karpeev /* determine global number of subdomains and the max number of local subdomains */ 332f746d493SDmitry Karpeev inwork[0] = inwork[1] = osm->n; 333b1a0a8a3SJed Brown ierr = MPI_Allreduce(inwork,outwork,1,MPIU_2INT,PetscMaxSum_Op,((PetscObject)pc)->comm);CHKERRQ(ierr); 334f746d493SDmitry Karpeev osm->nmax = outwork[0]; 335f746d493SDmitry Karpeev osm->N = outwork[1]; 336b1a0a8a3SJed Brown } 3376a4f0f73SDmitry Karpeev if (!osm->iis){ 33893c1ec32SDmitry Karpeev /* 3396a4f0f73SDmitry Karpeev The local number of subdomains was set in PCGASMSetTotalSubdomains() or PCGASMSetSubdomains(), 3406a4f0f73SDmitry Karpeev but the actual subdomains have not been supplied (in PCGASMSetSubdomains()). 3416a4f0f73SDmitry Karpeev We create the requisite number of inner subdomains on PETSC_COMM_SELF (for now). 34293c1ec32SDmitry Karpeev */ 3436a4f0f73SDmitry Karpeev ierr = PCGASMCreateLocalSubdomains(pc->pmat,osm->overlap,osm->n,&osm->iis,&osm->ois);CHKERRQ(ierr); 344b1a0a8a3SJed Brown } 3456a4f0f73SDmitry Karpeev 346b1a0a8a3SJed Brown ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 347b1a0a8a3SJed Brown flg = PETSC_FALSE; 348f746d493SDmitry Karpeev ierr = PetscOptionsGetBool(prefix,"-pc_gasm_print_subdomains",&flg,PETSC_NULL);CHKERRQ(ierr); 349f746d493SDmitry Karpeev if (flg) { ierr = PCGASMPrintSubdomains(pc);CHKERRQ(ierr); } 350b1a0a8a3SJed Brown 351b1a0a8a3SJed Brown if (osm->sort_indices) { 352f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 3536a4f0f73SDmitry Karpeev ierr = ISSort(osm->ois[i]);CHKERRQ(ierr); 3546a4f0f73SDmitry Karpeev ierr = ISSort(osm->iis[i]);CHKERRQ(ierr); 355b1a0a8a3SJed Brown } 356b1a0a8a3SJed Brown } 3576a4f0f73SDmitry Karpeev /* 3586a4f0f73SDmitry Karpeev Merge the ISs, create merged vectors and restrictions. 3596a4f0f73SDmitry Karpeev */ 3606a4f0f73SDmitry Karpeev /* Merge outer subdomain ISs and construct a restriction onto the disjoint union of local outer subdomains. */ 3616a4f0f73SDmitry Karpeev on = 0; 362f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 3636a4f0f73SDmitry Karpeev ierr = ISGetLocalSize(osm->ois[i],&oni);CHKERRQ(ierr); 3646a4f0f73SDmitry Karpeev on += oni; 365f746d493SDmitry Karpeev } 3666a4f0f73SDmitry Karpeev ierr = PetscMalloc(on*sizeof(PetscInt), &oidx);CHKERRQ(ierr); 3676a4f0f73SDmitry Karpeev on = 0; 368f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 3696a4f0f73SDmitry Karpeev ierr = ISGetLocalSize(osm->ois[i],&oni);CHKERRQ(ierr); 3706a4f0f73SDmitry Karpeev ierr = ISGetIndices(osm->ois[i],&oidxi);CHKERRQ(ierr); 3716a4f0f73SDmitry Karpeev ierr = PetscMemcpy(oidx+on, oidxi, sizeof(PetscInt)*oni);CHKERRQ(ierr); 3726a4f0f73SDmitry Karpeev ierr = ISRestoreIndices(osm->ois[i], &oidxi);CHKERRQ(ierr); 3736a4f0f73SDmitry Karpeev on += oni; 374f746d493SDmitry Karpeev } 3756a4f0f73SDmitry Karpeev ierr = ISCreateGeneral(((PetscObject)(pc))->comm, on, oidx, PETSC_OWN_POINTER, &gois);CHKERRQ(ierr); 376f746d493SDmitry Karpeev ierr = MatGetVecs(pc->pmat,&x,&y);CHKERRQ(ierr); 3776a4f0f73SDmitry Karpeev ierr = VecCreateMPI(((PetscObject)pc)->comm, on, PETSC_DECIDE, &osm->gx);CHKERRQ(ierr); 378f746d493SDmitry Karpeev ierr = VecDuplicate(osm->gx,&osm->gy);CHKERRQ(ierr); 3796a4f0f73SDmitry Karpeev ierr = VecGetOwnershipRange(osm->gx, &gofirst, PETSC_NULL);CHKERRQ(ierr); 3806a4f0f73SDmitry Karpeev ierr = ISCreateStride(((PetscObject)pc)->comm,on,gofirst,1, &goid);CHKERRQ(ierr); 3816a4f0f73SDmitry Karpeev ierr = VecScatterCreate(x,gois,osm->gx,goid, &(osm->gorestriction));CHKERRQ(ierr); 3826a4f0f73SDmitry Karpeev ierr = VecDestroy(&x); CHKERRQ(ierr); 383*7105d80bSDmitry Karpeev ierr = ISDestroy(&gois); CHKERRQ(ierr); 3846a4f0f73SDmitry Karpeev /* Merge inner subdomain ISs and construct a restriction onto the disjoint union of local inner subdomains. */ 3856a4f0f73SDmitry Karpeev { PetscInt ini; /* Number of indices the i-th a local inner subdomain. */ 3866a4f0f73SDmitry Karpeev PetscInt in; /* Number of indices in the disjoint uniont of local inner subdomains. */ 3876a4f0f73SDmitry Karpeev PetscInt *iidx; /* Global indices in the merged local inner subdomain. */ 3886a4f0f73SDmitry Karpeev PetscInt *ioidx; /* Global indices of the disjoint union of inner subdomains within the disjoint union of outer subdomains. */ 3896a4f0f73SDmitry Karpeev IS giis; /* IS for the disjoint union of inner subdomains. */ 3906a4f0f73SDmitry Karpeev IS giois; /* IS for the disjoint union of inner subdomains within the disjoint union of outer subdomains. */ 391f746d493SDmitry Karpeev /**/ 3926a4f0f73SDmitry Karpeev in = 0; 393f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 3946a4f0f73SDmitry Karpeev ierr = ISGetLocalSize(osm->iis[i],&ini);CHKERRQ(ierr); 3956a4f0f73SDmitry Karpeev in += ini; 396f746d493SDmitry Karpeev } 3976a4f0f73SDmitry Karpeev ierr = PetscMalloc(in*sizeof(PetscInt), &iidx); CHKERRQ(ierr); 3986a4f0f73SDmitry Karpeev ierr = PetscMalloc(in*sizeof(PetscInt), &ioidx);CHKERRQ(ierr); 3996a4f0f73SDmitry Karpeev ierr = VecGetOwnershipRange(osm->gx,&gofirst, PETSC_NULL); CHKERRQ(ierr); 4006a4f0f73SDmitry Karpeev in = 0; 4016a4f0f73SDmitry Karpeev on = 0; 402f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 4036a4f0f73SDmitry Karpeev const PetscInt *iidxi; /* Global indices of the i-th local inner subdomain. */ 4046a4f0f73SDmitry Karpeev ISLocalToGlobalMapping ltogi; /* Map from global to local indices of the i-th outer local subdomain. */ 4056a4f0f73SDmitry Karpeev PetscInt *ioidxi; /* Local indices of the i-th local inner subdomain within the local outer subdomain. */ 4066a4f0f73SDmitry Karpeev PetscInt ioni; /* Number of indices in ioidxi; if ioni != ini the inner subdomain is not a subdomain of the outer subdomain (error). */ 4076a4f0f73SDmitry Karpeev PetscInt k; 4086a4f0f73SDmitry Karpeev ierr = ISGetLocalSize(osm->iis[i],&ini);CHKERRQ(ierr); 4096a4f0f73SDmitry Karpeev ierr = ISGetLocalSize(osm->ois[i],&oni);CHKERRQ(ierr); 4106a4f0f73SDmitry Karpeev ierr = ISGetIndices(osm->iis[i],&iidxi);CHKERRQ(ierr); 4116a4f0f73SDmitry Karpeev ierr = PetscMemcpy(iidx+in, iidxi, sizeof(PetscInt)*ini);CHKERRQ(ierr); 4126a4f0f73SDmitry Karpeev ierr = ISLocalToGlobalMappingCreateIS(osm->ois[i],<ogi);CHKERRQ(ierr); 4136a4f0f73SDmitry Karpeev ioidxi = ioidx+in; 4146a4f0f73SDmitry Karpeev ierr = ISGlobalToLocalMappingApply(ltogi,IS_GTOLM_DROP,ini,iidxi,&ioni,ioidxi);CHKERRQ(ierr); 4156a4f0f73SDmitry Karpeev ierr = ISLocalToGlobalMappingDestroy(<ogi);CHKERRQ(ierr); 4166a4f0f73SDmitry Karpeev ierr = ISRestoreIndices(osm->iis[i], &iidxi);CHKERRQ(ierr); 4176a4f0f73SDmitry 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); 4186a4f0f73SDmitry Karpeev for(k = 0; k < ini; ++k) { 4196a4f0f73SDmitry Karpeev ioidxi[k] += gofirst+on; 420f746d493SDmitry Karpeev } 4216a4f0f73SDmitry Karpeev in += ini; 4226a4f0f73SDmitry Karpeev on += oni; 423f746d493SDmitry Karpeev } 4246a4f0f73SDmitry Karpeev ierr = ISCreateGeneral(((PetscObject)pc)->comm, in, iidx, PETSC_OWN_POINTER, &giis);CHKERRQ(ierr); 4256a4f0f73SDmitry Karpeev ierr = ISCreateGeneral(((PetscObject)pc)->comm, in, ioidx, PETSC_OWN_POINTER, &giois);CHKERRQ(ierr); 4266a4f0f73SDmitry Karpeev ierr = VecScatterCreate(y,giis,osm->gy,giois,&osm->girestriction);CHKERRQ(ierr); 4276a4f0f73SDmitry Karpeev ierr = VecDestroy(&y); CHKERRQ(ierr); 4286a4f0f73SDmitry Karpeev ierr = ISDestroy(&giis); CHKERRQ(ierr); 4296a4f0f73SDmitry Karpeev ierr = ISDestroy(&giois); CHKERRQ(ierr); 430b1a0a8a3SJed Brown } 4316a4f0f73SDmitry Karpeev ierr = ISDestroy(&goid);CHKERRQ(ierr); 432f746d493SDmitry Karpeev /* Create the subdomain work vectors. */ 433f746d493SDmitry Karpeev ierr = PetscMalloc(osm->n*sizeof(Vec),&osm->x);CHKERRQ(ierr); 434f746d493SDmitry Karpeev ierr = PetscMalloc(osm->n*sizeof(Vec),&osm->y);CHKERRQ(ierr); 435f746d493SDmitry Karpeev ierr = VecGetArray(osm->gx, &gxarray);CHKERRQ(ierr); 436f746d493SDmitry Karpeev ierr = VecGetArray(osm->gy, &gyarray);CHKERRQ(ierr); 4376a4f0f73SDmitry Karpeev for (i=0, on=0; i<osm->n; ++i, on += oni) { 4386a4f0f73SDmitry Karpeev PetscInt oNi; 4396a4f0f73SDmitry Karpeev ierr = ISGetLocalSize(osm->ois[i],&oni);CHKERRQ(ierr); 4406a4f0f73SDmitry Karpeev ierr = ISGetSize(osm->ois[i],&oNi);CHKERRQ(ierr); 4416a4f0f73SDmitry Karpeev ierr = VecCreateMPIWithArray(((PetscObject)(osm->ois[i]))->comm,1,oni,oNi,gxarray+on,&osm->x[i]);CHKERRQ(ierr); 4426a4f0f73SDmitry Karpeev ierr = VecCreateMPIWithArray(((PetscObject)(osm->ois[i]))->comm,1,oni,oNi,gyarray+on,&osm->y[i]);CHKERRQ(ierr); 443b1a0a8a3SJed Brown } 444f746d493SDmitry Karpeev ierr = VecRestoreArray(osm->gx, &gxarray);CHKERRQ(ierr); 445f746d493SDmitry Karpeev ierr = VecRestoreArray(osm->gy, &gyarray);CHKERRQ(ierr); 446b1a0a8a3SJed Brown /* Create the local solvers */ 447f746d493SDmitry Karpeev ierr = PetscMalloc(osm->n*sizeof(KSP *),&osm->ksp);CHKERRQ(ierr); 44844fe18e5SDmitry Karpeev for (i=0; i<osm->n; i++) { /* KSPs are local */ 4496a4f0f73SDmitry Karpeev ierr = KSPCreate(((PetscObject)(osm->ois[i]))->comm,&ksp);CHKERRQ(ierr); 450b1a0a8a3SJed Brown ierr = PetscLogObjectParent(pc,ksp);CHKERRQ(ierr); 451b1a0a8a3SJed Brown ierr = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,1);CHKERRQ(ierr); 452b1a0a8a3SJed Brown ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr); 453b1a0a8a3SJed Brown ierr = KSPGetPC(ksp,&subpc);CHKERRQ(ierr); 454b1a0a8a3SJed Brown ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 455b1a0a8a3SJed Brown ierr = KSPSetOptionsPrefix(ksp,prefix);CHKERRQ(ierr); 456b1a0a8a3SJed Brown ierr = KSPAppendOptionsPrefix(ksp,"sub_");CHKERRQ(ierr); 457b1a0a8a3SJed Brown osm->ksp[i] = ksp; 458b1a0a8a3SJed Brown } 459b1a0a8a3SJed Brown scall = MAT_INITIAL_MATRIX; 460b1a0a8a3SJed Brown 4616a4f0f73SDmitry Karpeev }/*if(!pc->setupcalled)*/ 4626a4f0f73SDmitry Karpeev else { 463b1a0a8a3SJed Brown /* 464b1a0a8a3SJed Brown Destroy the blocks from the previous iteration 465b1a0a8a3SJed Brown */ 466b1a0a8a3SJed Brown if (pc->flag == DIFFERENT_NONZERO_PATTERN) { 467f746d493SDmitry Karpeev ierr = MatDestroyMatrices(osm->n,&osm->pmat);CHKERRQ(ierr); 468b1a0a8a3SJed Brown scall = MAT_INITIAL_MATRIX; 469b1a0a8a3SJed Brown } 470b1a0a8a3SJed Brown } 471b1a0a8a3SJed Brown 472b1a0a8a3SJed Brown /* 473f746d493SDmitry Karpeev Extract out the submatrices. 474b1a0a8a3SJed Brown */ 47581a5c4aaSDmitry Karpeev if(size > 1) { 4766a4f0f73SDmitry Karpeev ierr = MatGetSubMatricesParallel(pc->pmat,osm->n,osm->ois, osm->ois,scall,&osm->pmat);CHKERRQ(ierr); 47781a5c4aaSDmitry Karpeev } 47881a5c4aaSDmitry Karpeev else { 4796a4f0f73SDmitry Karpeev ierr = MatGetSubMatrices(pc->pmat,osm->n,osm->ois, osm->ois,scall,&osm->pmat);CHKERRQ(ierr); 48081a5c4aaSDmitry Karpeev } 481b1a0a8a3SJed Brown if (scall == MAT_INITIAL_MATRIX) { 482b1a0a8a3SJed Brown ierr = PetscObjectGetOptionsPrefix((PetscObject)pc->pmat,&pprefix);CHKERRQ(ierr); 483f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 484b1a0a8a3SJed Brown ierr = PetscLogObjectParent(pc,osm->pmat[i]);CHKERRQ(ierr); 485b1a0a8a3SJed Brown ierr = PetscObjectSetOptionsPrefix((PetscObject)osm->pmat[i],pprefix);CHKERRQ(ierr); 486b1a0a8a3SJed Brown } 487b1a0a8a3SJed Brown } 488b1a0a8a3SJed Brown 489b1a0a8a3SJed Brown /* Return control to the user so that the submatrices can be modified (e.g., to apply 490b1a0a8a3SJed Brown different boundary conditions for the submatrices than for the global problem) */ 4916a4f0f73SDmitry Karpeev ierr = PCModifySubMatrices(pc,osm->n,osm->ois,osm->ois,osm->pmat,pc->modifysubmatricesP);CHKERRQ(ierr); 492b1a0a8a3SJed Brown 493b1a0a8a3SJed Brown /* 4946a4f0f73SDmitry Karpeev Loop over submatrices putting them into local ksps 495b1a0a8a3SJed Brown */ 496f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 497b1a0a8a3SJed Brown ierr = KSPSetOperators(osm->ksp[i],osm->pmat[i],osm->pmat[i],pc->flag);CHKERRQ(ierr); 498b1a0a8a3SJed Brown if (!pc->setupcalled) { 499b1a0a8a3SJed Brown ierr = KSPSetFromOptions(osm->ksp[i]);CHKERRQ(ierr); 500b1a0a8a3SJed Brown } 501b1a0a8a3SJed Brown } 502b1a0a8a3SJed Brown 503b1a0a8a3SJed Brown PetscFunctionReturn(0); 504b1a0a8a3SJed Brown } 505b1a0a8a3SJed Brown 506b1a0a8a3SJed Brown #undef __FUNCT__ 507f746d493SDmitry Karpeev #define __FUNCT__ "PCSetUpOnBlocks_GASM" 508f746d493SDmitry Karpeev static PetscErrorCode PCSetUpOnBlocks_GASM(PC pc) 509b1a0a8a3SJed Brown { 510f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 511b1a0a8a3SJed Brown PetscErrorCode ierr; 512b1a0a8a3SJed Brown PetscInt i; 513b1a0a8a3SJed Brown 514b1a0a8a3SJed Brown PetscFunctionBegin; 515f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 516b1a0a8a3SJed Brown ierr = KSPSetUp(osm->ksp[i]);CHKERRQ(ierr); 517b1a0a8a3SJed Brown } 518b1a0a8a3SJed Brown PetscFunctionReturn(0); 519b1a0a8a3SJed Brown } 520b1a0a8a3SJed Brown 521b1a0a8a3SJed Brown #undef __FUNCT__ 522f746d493SDmitry Karpeev #define __FUNCT__ "PCApply_GASM" 523f746d493SDmitry Karpeev static PetscErrorCode PCApply_GASM(PC pc,Vec x,Vec y) 524b1a0a8a3SJed Brown { 525f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 526b1a0a8a3SJed Brown PetscErrorCode ierr; 527f746d493SDmitry Karpeev PetscInt i; 528b1a0a8a3SJed Brown ScatterMode forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE; 529b1a0a8a3SJed Brown 530b1a0a8a3SJed Brown PetscFunctionBegin; 531b1a0a8a3SJed Brown /* 5326a4f0f73SDmitry Karpeev Support for limiting the restriction or interpolation only to the inner 533b1a0a8a3SJed Brown subdomain values (leaving the other values 0). 534b1a0a8a3SJed Brown */ 535f746d493SDmitry Karpeev if(!(osm->type & PC_GASM_RESTRICT)) { 536b1a0a8a3SJed Brown /* have to zero the work RHS since scatter may leave some slots empty */ 537f746d493SDmitry Karpeev ierr = VecZeroEntries(osm->gx);CHKERRQ(ierr); 5386a4f0f73SDmitry Karpeev ierr = VecScatterBegin(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 539b1a0a8a3SJed Brown } 5406a4f0f73SDmitry Karpeev else { 5416a4f0f73SDmitry Karpeev ierr = VecScatterBegin(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 542b1a0a8a3SJed Brown } 5436a4f0f73SDmitry Karpeev ierr = VecZeroEntries(osm->gy);CHKERRQ(ierr); 5446a4f0f73SDmitry Karpeev if(!(osm->type & PC_GASM_RESTRICT)) { 5456a4f0f73SDmitry Karpeev ierr = VecScatterEnd(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 5466a4f0f73SDmitry Karpeev } 5476a4f0f73SDmitry Karpeev else { 5486a4f0f73SDmitry Karpeev ierr = VecScatterEnd(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 5496a4f0f73SDmitry Karpeev } 55086b96d36SDmitry Karpeev /* do the subdomain solves */ 55186b96d36SDmitry Karpeev for (i=0; i<osm->n; ++i) { 552b1a0a8a3SJed Brown ierr = KSPSolve(osm->ksp[i],osm->x[i],osm->y[i]);CHKERRQ(ierr); 553b1a0a8a3SJed Brown } 5546a4f0f73SDmitry Karpeev ierr = VecZeroEntries(y); CHKERRQ(ierr); 5556a4f0f73SDmitry Karpeev if(!(osm->type & PC_GASM_INTERPOLATE)) { 5566a4f0f73SDmitry Karpeev ierr = VecScatterBegin(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); 5576a4f0f73SDmitry Karpeev ierr = VecScatterEnd(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); PetscFunctionReturn(0); 5586a4f0f73SDmitry Karpeev } 5596a4f0f73SDmitry Karpeev else { 5606a4f0f73SDmitry Karpeev ierr = VecScatterBegin(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); 5616a4f0f73SDmitry Karpeev ierr = VecScatterEnd(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); PetscFunctionReturn(0); 5626a4f0f73SDmitry Karpeev } 563b1a0a8a3SJed Brown } 564b1a0a8a3SJed Brown 565b1a0a8a3SJed Brown #undef __FUNCT__ 566f746d493SDmitry Karpeev #define __FUNCT__ "PCApplyTranspose_GASM" 567f746d493SDmitry Karpeev static PetscErrorCode PCApplyTranspose_GASM(PC pc,Vec x,Vec y) 568b1a0a8a3SJed Brown { 569f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 570b1a0a8a3SJed Brown PetscErrorCode ierr; 571f746d493SDmitry Karpeev PetscInt i; 572b1a0a8a3SJed Brown ScatterMode forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE; 573b1a0a8a3SJed Brown 574b1a0a8a3SJed Brown PetscFunctionBegin; 575b1a0a8a3SJed Brown /* 576b1a0a8a3SJed Brown Support for limiting the restriction or interpolation to only local 577b1a0a8a3SJed Brown subdomain values (leaving the other values 0). 578b1a0a8a3SJed Brown 579f746d493SDmitry Karpeev Note: these are reversed from the PCApply_GASM() because we are applying the 580b1a0a8a3SJed Brown transpose of the three terms 581b1a0a8a3SJed Brown */ 582f746d493SDmitry Karpeev if (!(osm->type & PC_GASM_INTERPOLATE)) { 583b1a0a8a3SJed Brown /* have to zero the work RHS since scatter may leave some slots empty */ 584f746d493SDmitry Karpeev ierr = VecZeroEntries(osm->gx);CHKERRQ(ierr); 5856a4f0f73SDmitry Karpeev ierr = VecScatterBegin(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 586b1a0a8a3SJed Brown } 5876a4f0f73SDmitry Karpeev else { 5886a4f0f73SDmitry Karpeev ierr = VecScatterBegin(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 589b1a0a8a3SJed Brown } 5906a4f0f73SDmitry Karpeev ierr = VecZeroEntries(osm->gy);CHKERRQ(ierr); 5916a4f0f73SDmitry Karpeev if (!(osm->type & PC_GASM_INTERPOLATE)) { 5926a4f0f73SDmitry Karpeev ierr = VecScatterEnd(osm->girestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 5936a4f0f73SDmitry Karpeev } 5946a4f0f73SDmitry Karpeev else { 5956a4f0f73SDmitry Karpeev ierr = VecScatterEnd(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 5966a4f0f73SDmitry Karpeev } 597b1a0a8a3SJed Brown /* do the local solves */ 598ab3e923fSDmitry 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. */ 599b1a0a8a3SJed Brown ierr = KSPSolveTranspose(osm->ksp[i],osm->x[i],osm->y[i]);CHKERRQ(ierr); 600b1a0a8a3SJed Brown } 6016a4f0f73SDmitry Karpeev ierr = VecZeroEntries(y); CHKERRQ(ierr); 6026a4f0f73SDmitry Karpeev if (!(osm->type & PC_GASM_RESTRICT)) { 6036a4f0f73SDmitry Karpeev ierr = VecScatterBegin(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); 6046a4f0f73SDmitry Karpeev ierr = VecScatterEnd(osm->girestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); 6056a4f0f73SDmitry Karpeev } 6066a4f0f73SDmitry Karpeev else { 6076a4f0f73SDmitry Karpeev ierr = VecScatterBegin(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); 6086a4f0f73SDmitry Karpeev ierr = VecScatterEnd(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); 6096a4f0f73SDmitry Karpeev } 6106a4f0f73SDmitry Karpeev 611b1a0a8a3SJed Brown PetscFunctionReturn(0); 612b1a0a8a3SJed Brown } 613b1a0a8a3SJed Brown 614b1a0a8a3SJed Brown #undef __FUNCT__ 615a06653b4SBarry Smith #define __FUNCT__ "PCReset_GASM" 616a06653b4SBarry Smith static PetscErrorCode PCReset_GASM(PC pc) 617b1a0a8a3SJed Brown { 618f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 619b1a0a8a3SJed Brown PetscErrorCode ierr; 620b1a0a8a3SJed Brown PetscInt i; 621b1a0a8a3SJed Brown 622b1a0a8a3SJed Brown PetscFunctionBegin; 623b1a0a8a3SJed Brown if (osm->ksp) { 624f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 625a06653b4SBarry Smith ierr = KSPReset(osm->ksp[i]);CHKERRQ(ierr); 626b1a0a8a3SJed Brown } 627b1a0a8a3SJed Brown } 628b1a0a8a3SJed Brown if (osm->pmat) { 629f746d493SDmitry Karpeev if (osm->n > 0) { 630ab3e923fSDmitry Karpeev ierr = MatDestroyMatrices(osm->n,&osm->pmat);CHKERRQ(ierr); 631b1a0a8a3SJed Brown } 632b1a0a8a3SJed Brown } 633d34fcf5fSBarry Smith if (osm->x) { 634f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 635fcfd50ebSBarry Smith ierr = VecDestroy(&osm->x[i]);CHKERRQ(ierr); 636fcfd50ebSBarry Smith ierr = VecDestroy(&osm->y[i]);CHKERRQ(ierr); 637b1a0a8a3SJed Brown } 638d34fcf5fSBarry Smith } 6396bf464f9SBarry Smith ierr = VecDestroy(&osm->gx);CHKERRQ(ierr); 6406bf464f9SBarry Smith ierr = VecDestroy(&osm->gy);CHKERRQ(ierr); 641ab3e923fSDmitry Karpeev 6426a4f0f73SDmitry Karpeev ierr = VecScatterDestroy(&osm->gorestriction);CHKERRQ(ierr); 6436a4f0f73SDmitry Karpeev ierr = VecScatterDestroy(&osm->girestriction);CHKERRQ(ierr); 6446a4f0f73SDmitry Karpeev if (osm->ois) { 6456a4f0f73SDmitry Karpeev ierr = PCGASMDestroySubdomains(osm->n,osm->ois,osm->iis);CHKERRQ(ierr); 6466a4f0f73SDmitry Karpeev osm->ois = 0; 6476a4f0f73SDmitry Karpeev osm->iis = 0; 648b6b0974bSDmitry Karpeev } 649a06653b4SBarry Smith PetscFunctionReturn(0); 650a06653b4SBarry Smith } 651a06653b4SBarry Smith 652a06653b4SBarry Smith #undef __FUNCT__ 653a06653b4SBarry Smith #define __FUNCT__ "PCDestroy_GASM" 654a06653b4SBarry Smith static PetscErrorCode PCDestroy_GASM(PC pc) 655a06653b4SBarry Smith { 656a06653b4SBarry Smith PC_GASM *osm = (PC_GASM*)pc->data; 657a06653b4SBarry Smith PetscErrorCode ierr; 658a06653b4SBarry Smith PetscInt i; 659a06653b4SBarry Smith 660a06653b4SBarry Smith PetscFunctionBegin; 661a06653b4SBarry Smith ierr = PCReset_GASM(pc);CHKERRQ(ierr); 662a06653b4SBarry Smith if (osm->ksp) { 663a06653b4SBarry Smith for (i=0; i<osm->n; i++) { 6646bf464f9SBarry Smith ierr = KSPDestroy(&osm->ksp[i]);CHKERRQ(ierr); 665a06653b4SBarry Smith } 666a06653b4SBarry Smith ierr = PetscFree(osm->ksp);CHKERRQ(ierr); 667a06653b4SBarry Smith } 668a06653b4SBarry Smith ierr = PetscFree(osm->x);CHKERRQ(ierr); 669a06653b4SBarry Smith ierr = PetscFree(osm->y);CHKERRQ(ierr); 670c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 671b1a0a8a3SJed Brown PetscFunctionReturn(0); 672b1a0a8a3SJed Brown } 673b1a0a8a3SJed Brown 674b1a0a8a3SJed Brown #undef __FUNCT__ 675f746d493SDmitry Karpeev #define __FUNCT__ "PCSetFromOptions_GASM" 676ab3e923fSDmitry Karpeev static PetscErrorCode PCSetFromOptions_GASM(PC pc) { 677f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 678b1a0a8a3SJed Brown PetscErrorCode ierr; 679b1a0a8a3SJed Brown PetscInt blocks,ovl; 680b1a0a8a3SJed Brown PetscBool symset,flg; 681f746d493SDmitry Karpeev PCGASMType gasmtype; 682b1a0a8a3SJed Brown 683b1a0a8a3SJed Brown PetscFunctionBegin; 684b1a0a8a3SJed Brown /* set the type to symmetric if matrix is symmetric */ 685b1a0a8a3SJed Brown if (!osm->type_set && pc->pmat) { 686b1a0a8a3SJed Brown ierr = MatIsSymmetricKnown(pc->pmat,&symset,&flg);CHKERRQ(ierr); 687f746d493SDmitry Karpeev if (symset && flg) { osm->type = PC_GASM_BASIC; } 688b1a0a8a3SJed Brown } 68944fe18e5SDmitry Karpeev ierr = PetscOptionsHead("Generalized additive Schwarz options");CHKERRQ(ierr); 6906a4f0f73SDmitry Karpeev ierr = PetscOptionsInt("-pc_gasm_total_subdomains","Total number of subdomains across communicator","PCGASMSetTotalSubdomains",osm->n,&blocks,&flg);CHKERRQ(ierr); 6916a4f0f73SDmitry Karpeev osm->create_local = PETSC_TRUE; 6926a4f0f73SDmitry 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); 6936a4f0f73SDmitry Karpeev if(!osm->create_local) SETERRQ(((PetscObject)pc)->comm, PETSC_ERR_SUP, "No support for autocreation of nonlocal subdomains yet."); 6946a4f0f73SDmitry Karpeev 6956a4f0f73SDmitry Karpeev if (flg) {ierr = PCGASMSetTotalSubdomains(pc,blocks,osm->create_local);CHKERRQ(ierr); } 69606b43e7eSDmitry Karpeev ierr = PetscOptionsInt("-pc_gasm_overlap","Number of overlapping degrees of freedom","PCGASMSetOverlap",osm->overlap,&ovl,&flg);CHKERRQ(ierr); 697f746d493SDmitry Karpeev if (flg) {ierr = PCGASMSetOverlap(pc,ovl);CHKERRQ(ierr); } 698b1a0a8a3SJed Brown flg = PETSC_FALSE; 699f746d493SDmitry Karpeev ierr = PetscOptionsEnum("-pc_gasm_type","Type of restriction/extension","PCGASMSetType",PCGASMTypes,(PetscEnum)osm->type,(PetscEnum*)&gasmtype,&flg);CHKERRQ(ierr); 700f746d493SDmitry Karpeev if (flg) {ierr = PCGASMSetType(pc,gasmtype);CHKERRQ(ierr); } 701b1a0a8a3SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 702b1a0a8a3SJed Brown PetscFunctionReturn(0); 703b1a0a8a3SJed Brown } 704b1a0a8a3SJed Brown 705b1a0a8a3SJed Brown /*------------------------------------------------------------------------------------*/ 706b1a0a8a3SJed Brown 707b1a0a8a3SJed Brown EXTERN_C_BEGIN 708b1a0a8a3SJed Brown #undef __FUNCT__ 70906b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMSetSubdomains_GASM" 7106a4f0f73SDmitry Karpeev PetscErrorCode PCGASMSetSubdomains_GASM(PC pc,PetscInt n,IS is_local[],IS is[]) 711b1a0a8a3SJed Brown { 712f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 713b1a0a8a3SJed Brown PetscErrorCode ierr; 714b1a0a8a3SJed Brown PetscInt i; 715b1a0a8a3SJed Brown 716b1a0a8a3SJed Brown PetscFunctionBegin; 71769ca1f37SDmitry Karpeev if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Each process must have 1 or more subdomains, n = %D",n); 7186a4f0f73SDmitry Karpeev if (pc->setupcalled && (n != osm->n || is_local)) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetSubdomains() should be called before calling PCSetUp()."); 719b1a0a8a3SJed Brown 720b1a0a8a3SJed Brown if (!pc->setupcalled) { 72193c1ec32SDmitry Karpeev osm->n = n; 7226a4f0f73SDmitry Karpeev osm->ois = 0; 7236a4f0f73SDmitry Karpeev osm->iis = 0; 724b1a0a8a3SJed Brown if (is) { 725b1a0a8a3SJed Brown for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)is[i]);CHKERRQ(ierr);} 726b1a0a8a3SJed Brown } 727b1a0a8a3SJed Brown if (is_local) { 728b1a0a8a3SJed Brown for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)is_local[i]);CHKERRQ(ierr);} 729b1a0a8a3SJed Brown } 7306a4f0f73SDmitry Karpeev if (osm->ois) { 7316a4f0f73SDmitry Karpeev ierr = PCGASMDestroySubdomains(osm->n,osm->iis,osm->ois);CHKERRQ(ierr); 732b1a0a8a3SJed Brown } 733b1a0a8a3SJed Brown if (is) { 7346a4f0f73SDmitry Karpeev ierr = PetscMalloc(n*sizeof(IS),&osm->ois);CHKERRQ(ierr); 7356a4f0f73SDmitry Karpeev for (i=0; i<n; i++) { osm->ois[i] = is[i]; } 7366a4f0f73SDmitry Karpeev /* Flag indicating that the user has set outer subdomains, so PCGASM should not increase their size. */ 737b1a0a8a3SJed Brown osm->overlap = -1; 738b1a0a8a3SJed Brown } 739b1a0a8a3SJed Brown if (is_local) { 7406a4f0f73SDmitry Karpeev ierr = PetscMalloc(n*sizeof(IS),&osm->iis);CHKERRQ(ierr); 7416a4f0f73SDmitry Karpeev for (i=0; i<n; i++) { osm->iis[i] = is_local[i]; } 742b1a0a8a3SJed Brown } 743b1a0a8a3SJed Brown } 744b1a0a8a3SJed Brown PetscFunctionReturn(0); 745b1a0a8a3SJed Brown } 746b1a0a8a3SJed Brown EXTERN_C_END 747b1a0a8a3SJed Brown 748b1a0a8a3SJed Brown EXTERN_C_BEGIN 749b1a0a8a3SJed Brown #undef __FUNCT__ 7506a4f0f73SDmitry Karpeev #define __FUNCT__ "PCGASMSetTotalSubdomains_GASM" 7516a4f0f73SDmitry Karpeev PetscErrorCode PCGASMSetTotalSubdomains_GASM(PC pc,PetscInt N, PetscBool create_local) { 752f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 753b1a0a8a3SJed Brown PetscErrorCode ierr; 754b1a0a8a3SJed Brown PetscMPIInt rank,size; 755b1a0a8a3SJed Brown PetscInt n; 7566a4f0f73SDmitry Karpeev PetscInt Nmin, Nmax; 757b1a0a8a3SJed Brown PetscFunctionBegin; 7586a4f0f73SDmitry Karpeev if(!create_local) SETERRQ(((PetscObject)pc)->comm, PETSC_ERR_SUP, "No suppor for autocreation of nonlocal subdomains."); 7596a4f0f73SDmitry Karpeev if (N < 1) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Total number of subdomains must be > 0, N = %D",N); 7606a4f0f73SDmitry Karpeev ierr = MPI_Allreduce(&N,&Nmin,1,MPI_INT,MPI_MIN,((PetscObject)pc)->comm); CHKERRQ(ierr); 7616a4f0f73SDmitry Karpeev ierr = MPI_Allreduce(&N,&Nmax,1,MPI_INT,MPI_MAX,((PetscObject)pc)->comm); CHKERRQ(ierr); 7626a4f0f73SDmitry Karpeev if(Nmin != Nmax) 7636a4f0f73SDmitry 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); 764b1a0a8a3SJed Brown 7656a4f0f73SDmitry Karpeev osm->create_local = create_local; 766b1a0a8a3SJed Brown /* 767b1a0a8a3SJed Brown Split the subdomains equally among all processors 768b1a0a8a3SJed Brown */ 769b1a0a8a3SJed Brown ierr = MPI_Comm_rank(((PetscObject)pc)->comm,&rank);CHKERRQ(ierr); 770b1a0a8a3SJed Brown ierr = MPI_Comm_size(((PetscObject)pc)->comm,&size);CHKERRQ(ierr); 771b1a0a8a3SJed Brown n = N/size + ((N % size) > rank); 7726a4f0f73SDmitry 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); 773f746d493SDmitry Karpeev if (pc->setupcalled && n != osm->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetTotalSubdomains() should be called before PCSetUp()."); 774b1a0a8a3SJed Brown if (!pc->setupcalled) { 7756a4f0f73SDmitry Karpeev if (osm->ois) { 7766a4f0f73SDmitry Karpeev ierr = PCGASMDestroySubdomains(osm->n,osm->iis,osm->ois);CHKERRQ(ierr); 777b1a0a8a3SJed Brown } 77893c1ec32SDmitry Karpeev osm->N = N; 779f746d493SDmitry Karpeev osm->n = n; 78006b43e7eSDmitry Karpeev osm->nmax = N/size + ((N%size)?1:0); 7816a4f0f73SDmitry Karpeev osm->ois = 0; 7826a4f0f73SDmitry Karpeev osm->iis = 0; 783b1a0a8a3SJed Brown } 784b1a0a8a3SJed Brown PetscFunctionReturn(0); 78506b43e7eSDmitry Karpeev } 786b1a0a8a3SJed Brown EXTERN_C_END 787b1a0a8a3SJed Brown 788b1a0a8a3SJed Brown EXTERN_C_BEGIN 789b1a0a8a3SJed Brown #undef __FUNCT__ 790f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetOverlap_GASM" 7917087cfbeSBarry Smith PetscErrorCode PCGASMSetOverlap_GASM(PC pc,PetscInt ovl) 792b1a0a8a3SJed Brown { 793f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 794b1a0a8a3SJed Brown 795b1a0a8a3SJed Brown PetscFunctionBegin; 796b1a0a8a3SJed Brown if (ovl < 0) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap value requested"); 797f746d493SDmitry Karpeev if (pc->setupcalled && ovl != osm->overlap) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetOverlap() should be called before PCSetUp()."); 798b1a0a8a3SJed Brown if (!pc->setupcalled) { 799b1a0a8a3SJed Brown osm->overlap = ovl; 800b1a0a8a3SJed Brown } 801b1a0a8a3SJed Brown PetscFunctionReturn(0); 802b1a0a8a3SJed Brown } 803b1a0a8a3SJed Brown EXTERN_C_END 804b1a0a8a3SJed Brown 805b1a0a8a3SJed Brown EXTERN_C_BEGIN 806b1a0a8a3SJed Brown #undef __FUNCT__ 807f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetType_GASM" 8087087cfbeSBarry Smith PetscErrorCode PCGASMSetType_GASM(PC pc,PCGASMType type) 809b1a0a8a3SJed Brown { 810f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 811b1a0a8a3SJed Brown 812b1a0a8a3SJed Brown PetscFunctionBegin; 813b1a0a8a3SJed Brown osm->type = type; 814b1a0a8a3SJed Brown osm->type_set = PETSC_TRUE; 815b1a0a8a3SJed Brown PetscFunctionReturn(0); 816b1a0a8a3SJed Brown } 817b1a0a8a3SJed Brown EXTERN_C_END 818b1a0a8a3SJed Brown 819b1a0a8a3SJed Brown EXTERN_C_BEGIN 820b1a0a8a3SJed Brown #undef __FUNCT__ 821f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetSortIndices_GASM" 8227087cfbeSBarry Smith PetscErrorCode PCGASMSetSortIndices_GASM(PC pc,PetscBool doSort) 823b1a0a8a3SJed Brown { 824f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 825b1a0a8a3SJed Brown 826b1a0a8a3SJed Brown PetscFunctionBegin; 827b1a0a8a3SJed Brown osm->sort_indices = doSort; 828b1a0a8a3SJed Brown PetscFunctionReturn(0); 829b1a0a8a3SJed Brown } 830b1a0a8a3SJed Brown EXTERN_C_END 831b1a0a8a3SJed Brown 832b1a0a8a3SJed Brown EXTERN_C_BEGIN 833b1a0a8a3SJed Brown #undef __FUNCT__ 834f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMGetSubKSP_GASM" 83544fe18e5SDmitry Karpeev /* 83644fe18e5SDmitry Karpeev FIX: This routine might need to be modified once multiple ranks per subdomain are allowed. 83744fe18e5SDmitry Karpeev In particular, it would upset the global subdomain number calculation. 83844fe18e5SDmitry Karpeev */ 8397087cfbeSBarry Smith PetscErrorCode PCGASMGetSubKSP_GASM(PC pc,PetscInt *n,PetscInt *first,KSP **ksp) 840b1a0a8a3SJed Brown { 841f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 842b1a0a8a3SJed Brown PetscErrorCode ierr; 843b1a0a8a3SJed Brown 844b1a0a8a3SJed Brown PetscFunctionBegin; 845ab3e923fSDmitry 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"); 846b1a0a8a3SJed Brown 847ab3e923fSDmitry Karpeev if (n) { 848ab3e923fSDmitry Karpeev *n = osm->n; 849b1a0a8a3SJed Brown } 850ab3e923fSDmitry Karpeev if (first) { 851ab3e923fSDmitry Karpeev ierr = MPI_Scan(&osm->n,first,1,MPIU_INT,MPI_SUM,((PetscObject)pc)->comm);CHKERRQ(ierr); 852ab3e923fSDmitry Karpeev *first -= osm->n; 853b1a0a8a3SJed Brown } 854b1a0a8a3SJed Brown if (ksp) { 855b1a0a8a3SJed Brown /* Assume that local solves are now different; not necessarily 85606b43e7eSDmitry Karpeev true, though! This flag is used only for PCView_GASM() */ 857b1a0a8a3SJed Brown *ksp = osm->ksp; 8586a4f0f73SDmitry Karpeev osm->same_subdomain_solvers = PETSC_FALSE; 859b1a0a8a3SJed Brown } 860b1a0a8a3SJed Brown PetscFunctionReturn(0); 861ab3e923fSDmitry Karpeev }/* PCGASMGetSubKSP_GASM() */ 862b1a0a8a3SJed Brown EXTERN_C_END 863b1a0a8a3SJed Brown 864b1a0a8a3SJed Brown 865b1a0a8a3SJed Brown #undef __FUNCT__ 86606b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMSetSubdomains" 867b1a0a8a3SJed Brown /*@C 86806b43e7eSDmitry Karpeev PCGASMSetSubdomains - Sets the subdomains for this processor 86906b43e7eSDmitry Karpeev for the additive Schwarz preconditioner. 870b1a0a8a3SJed Brown 871b1a0a8a3SJed Brown Collective on PC 872b1a0a8a3SJed Brown 873b1a0a8a3SJed Brown Input Parameters: 874b1a0a8a3SJed Brown + pc - the preconditioner context 87506b43e7eSDmitry Karpeev . n - the number of subdomains for this processor 8766a4f0f73SDmitry Karpeev . iis - the index sets that define this processor's local inner subdomains 877b1a0a8a3SJed Brown (or PETSC_NULL for PETSc to determine subdomains) 8786a4f0f73SDmitry Karpeev - ois- the index sets that define this processor's local outer subdomains 8796a4f0f73SDmitry Karpeev (or PETSC_NULL to use the same as iis) 880b1a0a8a3SJed Brown 881b1a0a8a3SJed Brown Notes: 88206b43e7eSDmitry Karpeev The IS indices use the parallel, global numbering of the vector entries. 8836a4f0f73SDmitry Karpeev Inner subdomains are those where the correction is applied. 8846a4f0f73SDmitry Karpeev Outer subdomains are those where the residual necessary to obtain the 8856a4f0f73SDmitry Karpeev corrections is obtained (see PCGASMType for the use of inner/outer subdomains). 8866a4f0f73SDmitry Karpeev Both inner and outer subdomains can extend over several processors. 8876a4f0f73SDmitry Karpeev This processor's portion of a subdomain is known as a local subdomain. 8886a4f0f73SDmitry Karpeev 8896a4f0f73SDmitry Karpeev By default the GASM preconditioner uses 1 (local) subdomain per processor. 8906a4f0f73SDmitry Karpeev Use PCGASMSetTotalSubdomains() to set the total number of subdomains across 8916a4f0f73SDmitry Karpeev all processors that PCGASM will create automatically, and to specify whether 8926a4f0f73SDmitry Karpeev they should be local or not. 8936a4f0f73SDmitry Karpeev 894b1a0a8a3SJed Brown 895b1a0a8a3SJed Brown Level: advanced 896b1a0a8a3SJed Brown 89706b43e7eSDmitry Karpeev .keywords: PC, GASM, set, subdomains, additive Schwarz 898b1a0a8a3SJed Brown 8996a4f0f73SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetOverlap(), PCGASMGetSubKSP(), 90006b43e7eSDmitry Karpeev PCGASMCreateSubdomains2D(), PCGASMGetSubdomains() 901b1a0a8a3SJed Brown @*/ 9026a4f0f73SDmitry Karpeev PetscErrorCode PCGASMSetSubdomains(PC pc,PetscInt n,IS iis[],IS ois[]) 903b1a0a8a3SJed Brown { 904b1a0a8a3SJed Brown PetscErrorCode ierr; 905b1a0a8a3SJed Brown 906b1a0a8a3SJed Brown PetscFunctionBegin; 907b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 9086a4f0f73SDmitry Karpeev ierr = PetscTryMethod(pc,"PCGASMSetSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,n,iis,ois));CHKERRQ(ierr); 909b1a0a8a3SJed Brown PetscFunctionReturn(0); 910b1a0a8a3SJed Brown } 911b1a0a8a3SJed Brown 912b1a0a8a3SJed Brown #undef __FUNCT__ 9136a4f0f73SDmitry Karpeev #define __FUNCT__ "PCGASMSetTotalSubdomains" 914b1a0a8a3SJed Brown /*@C 9156a4f0f73SDmitry Karpeev PCGASMSetTotalSubdomains - Sets the total number of subdomains to use in the generalized additive 9166a4f0f73SDmitry Karpeev Schwarz preconditioner. The number of subdomains is cumulative across all processors in pc's 9176a4f0f73SDmitry Karpeev communicator. Either all or no processors in the PC communicator must call this routine with 9186a4f0f73SDmitry Karpeev the same N. The subdomains will be created automatically during PCSetUp(). 919b1a0a8a3SJed Brown 920b1a0a8a3SJed Brown Collective on PC 921b1a0a8a3SJed Brown 922b1a0a8a3SJed Brown Input Parameters: 923b1a0a8a3SJed Brown + pc - the preconditioner context 9246a4f0f73SDmitry Karpeev . N - the total number of subdomains cumulative across all processors 9256a4f0f73SDmitry Karpeev - create_local - whether the subdomains to be created are to be local 926b1a0a8a3SJed Brown 927b1a0a8a3SJed Brown Options Database Key: 9286a4f0f73SDmitry Karpeev To set the total number of subdomains and let PCGASM autocreate them, rather than specify the index sets, use the following options: 9296a4f0f73SDmitry Karpeev + -pc_gasm_total_subdomains <n> - sets the total number of subdomains to be autocreated by PCGASM 9306a4f0f73SDmitry Karpeev - -pc_gasm_subdomains_create_local <true|false> - whether autocreated subdomains should be local or not (default is true) 931b1a0a8a3SJed Brown 93206b43e7eSDmitry Karpeev By default the GASM preconditioner uses 1 subdomain per processor. 933b1a0a8a3SJed Brown 9346a4f0f73SDmitry Karpeev 9356a4f0f73SDmitry Karpeev Use PCGASMSetSubdomains() to set subdomains explicitly or to set different numbers 9366a4f0f73SDmitry Karpeev of subdomains per processor. 937b1a0a8a3SJed Brown 938b1a0a8a3SJed Brown Level: advanced 939b1a0a8a3SJed Brown 940f746d493SDmitry Karpeev .keywords: PC, GASM, set, total, global, subdomains, additive Schwarz 941b1a0a8a3SJed Brown 94206b43e7eSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMSetOverlap(), PCGASMGetSubKSP(), 943f746d493SDmitry Karpeev PCGASMCreateSubdomains2D() 944b1a0a8a3SJed Brown @*/ 9456a4f0f73SDmitry Karpeev PetscErrorCode PCGASMSetTotalSubdomains(PC pc,PetscInt N, PetscBool create_local) 946b1a0a8a3SJed Brown { 947b1a0a8a3SJed Brown PetscErrorCode ierr; 948b1a0a8a3SJed Brown 949b1a0a8a3SJed Brown PetscFunctionBegin; 950b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 9516a4f0f73SDmitry Karpeev ierr = PetscTryMethod(pc,"PCGASMSetTotalSubdomains_C",(PC,PetscInt,PetscBool),(pc,N,create_local));CHKERRQ(ierr); 952b1a0a8a3SJed Brown PetscFunctionReturn(0); 953b1a0a8a3SJed Brown } 954b1a0a8a3SJed Brown 955b1a0a8a3SJed Brown #undef __FUNCT__ 956f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetOverlap" 957b1a0a8a3SJed Brown /*@ 958f746d493SDmitry Karpeev PCGASMSetOverlap - Sets the overlap between a pair of subdomains for the 959b1a0a8a3SJed Brown additive Schwarz preconditioner. Either all or no processors in the 960b1a0a8a3SJed Brown PC communicator must call this routine. 961b1a0a8a3SJed Brown 962b1a0a8a3SJed Brown Logically Collective on PC 963b1a0a8a3SJed Brown 964b1a0a8a3SJed Brown Input Parameters: 965b1a0a8a3SJed Brown + pc - the preconditioner context 966b1a0a8a3SJed Brown - ovl - the amount of overlap between subdomains (ovl >= 0, default value = 1) 967b1a0a8a3SJed Brown 968b1a0a8a3SJed Brown Options Database Key: 96906b43e7eSDmitry Karpeev . -pc_gasm_overlap <overlap> - Sets overlap 970b1a0a8a3SJed Brown 971b1a0a8a3SJed Brown Notes: 97206b43e7eSDmitry Karpeev By default the GASM preconditioner uses 1 subdomain per processor. To use 9736a4f0f73SDmitry Karpeev multiple subdomain per perocessor, see PCGASMSetTotalSubdomains() or 9746a4f0f73SDmitry Karpeev PCGASMSetSubdomains() (and the option -pc_gasm_total_subdomains <n>). 975b1a0a8a3SJed Brown 976b1a0a8a3SJed Brown The overlap defaults to 1, so if one desires that no additional 977b1a0a8a3SJed Brown overlap be computed beyond what may have been set with a call to 9786a4f0f73SDmitry Karpeev PCGASMSetTotalSubdomains() or PCGASMSetSubdomains(), then ovl 979b1a0a8a3SJed Brown must be set to be 0. In particular, if one does not explicitly set 98006b43e7eSDmitry Karpeev the subdomains in application code, then all overlap would be computed 981f746d493SDmitry Karpeev internally by PETSc, and using an overlap of 0 would result in an GASM 982b1a0a8a3SJed Brown variant that is equivalent to the block Jacobi preconditioner. 983b1a0a8a3SJed Brown 984b1a0a8a3SJed Brown Note that one can define initial index sets with any overlap via 98506b43e7eSDmitry Karpeev PCGASMSetSubdomains(); the routine PCGASMSetOverlap() merely allows 98606b43e7eSDmitry Karpeev PETSc to extend that overlap further, if desired. 987b1a0a8a3SJed Brown 988b1a0a8a3SJed Brown Level: intermediate 989b1a0a8a3SJed Brown 990f746d493SDmitry Karpeev .keywords: PC, GASM, set, overlap 991b1a0a8a3SJed Brown 9926a4f0f73SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetSubdomains(), PCGASMGetSubKSP(), 99306b43e7eSDmitry Karpeev PCGASMCreateSubdomains2D(), PCGASMGetSubdomains() 994b1a0a8a3SJed Brown @*/ 9957087cfbeSBarry Smith PetscErrorCode PCGASMSetOverlap(PC pc,PetscInt ovl) 996b1a0a8a3SJed Brown { 997b1a0a8a3SJed Brown PetscErrorCode ierr; 998b1a0a8a3SJed Brown 999b1a0a8a3SJed Brown PetscFunctionBegin; 1000b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1001b1a0a8a3SJed Brown PetscValidLogicalCollectiveInt(pc,ovl,2); 1002f746d493SDmitry Karpeev ierr = PetscTryMethod(pc,"PCGASMSetOverlap_C",(PC,PetscInt),(pc,ovl));CHKERRQ(ierr); 1003b1a0a8a3SJed Brown PetscFunctionReturn(0); 1004b1a0a8a3SJed Brown } 1005b1a0a8a3SJed Brown 1006b1a0a8a3SJed Brown #undef __FUNCT__ 1007f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetType" 1008b1a0a8a3SJed Brown /*@ 1009f746d493SDmitry Karpeev PCGASMSetType - Sets the type of restriction and interpolation used 1010b1a0a8a3SJed Brown for local problems in the additive Schwarz method. 1011b1a0a8a3SJed Brown 1012b1a0a8a3SJed Brown Logically Collective on PC 1013b1a0a8a3SJed Brown 1014b1a0a8a3SJed Brown Input Parameters: 1015b1a0a8a3SJed Brown + pc - the preconditioner context 1016f746d493SDmitry Karpeev - type - variant of GASM, one of 1017b1a0a8a3SJed Brown .vb 1018f746d493SDmitry Karpeev PC_GASM_BASIC - full interpolation and restriction 1019f746d493SDmitry Karpeev PC_GASM_RESTRICT - full restriction, local processor interpolation 1020f746d493SDmitry Karpeev PC_GASM_INTERPOLATE - full interpolation, local processor restriction 1021f746d493SDmitry Karpeev PC_GASM_NONE - local processor restriction and interpolation 1022b1a0a8a3SJed Brown .ve 1023b1a0a8a3SJed Brown 1024b1a0a8a3SJed Brown Options Database Key: 1025f746d493SDmitry Karpeev . -pc_gasm_type [basic,restrict,interpolate,none] - Sets GASM type 1026b1a0a8a3SJed Brown 1027b1a0a8a3SJed Brown Level: intermediate 1028b1a0a8a3SJed Brown 1029f746d493SDmitry Karpeev .keywords: PC, GASM, set, type 1030b1a0a8a3SJed Brown 10316a4f0f73SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetSubdomains(), PCGASMGetSubKSP(), 1032f746d493SDmitry Karpeev PCGASMCreateSubdomains2D() 1033b1a0a8a3SJed Brown @*/ 10347087cfbeSBarry Smith PetscErrorCode PCGASMSetType(PC pc,PCGASMType type) 1035b1a0a8a3SJed Brown { 1036b1a0a8a3SJed Brown PetscErrorCode ierr; 1037b1a0a8a3SJed Brown 1038b1a0a8a3SJed Brown PetscFunctionBegin; 1039b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1040b1a0a8a3SJed Brown PetscValidLogicalCollectiveEnum(pc,type,2); 1041f746d493SDmitry Karpeev ierr = PetscTryMethod(pc,"PCGASMSetType_C",(PC,PCGASMType),(pc,type));CHKERRQ(ierr); 1042b1a0a8a3SJed Brown PetscFunctionReturn(0); 1043b1a0a8a3SJed Brown } 1044b1a0a8a3SJed Brown 1045b1a0a8a3SJed Brown #undef __FUNCT__ 1046f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetSortIndices" 1047b1a0a8a3SJed Brown /*@ 1048f746d493SDmitry Karpeev PCGASMSetSortIndices - Determines whether subdomain indices are sorted. 1049b1a0a8a3SJed Brown 1050b1a0a8a3SJed Brown Logically Collective on PC 1051b1a0a8a3SJed Brown 1052b1a0a8a3SJed Brown Input Parameters: 1053b1a0a8a3SJed Brown + pc - the preconditioner context 1054b1a0a8a3SJed Brown - doSort - sort the subdomain indices 1055b1a0a8a3SJed Brown 1056b1a0a8a3SJed Brown Level: intermediate 1057b1a0a8a3SJed Brown 1058f746d493SDmitry Karpeev .keywords: PC, GASM, set, type 1059b1a0a8a3SJed Brown 10606a4f0f73SDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMSetTotalSubdomains(), PCGASMGetSubKSP(), 1061f746d493SDmitry Karpeev PCGASMCreateSubdomains2D() 1062b1a0a8a3SJed Brown @*/ 10637087cfbeSBarry Smith PetscErrorCode PCGASMSetSortIndices(PC pc,PetscBool doSort) 1064b1a0a8a3SJed Brown { 1065b1a0a8a3SJed Brown PetscErrorCode ierr; 1066b1a0a8a3SJed Brown 1067b1a0a8a3SJed Brown PetscFunctionBegin; 1068b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1069b1a0a8a3SJed Brown PetscValidLogicalCollectiveBool(pc,doSort,2); 1070f746d493SDmitry Karpeev ierr = PetscTryMethod(pc,"PCGASMSetSortIndices_C",(PC,PetscBool),(pc,doSort));CHKERRQ(ierr); 1071b1a0a8a3SJed Brown PetscFunctionReturn(0); 1072b1a0a8a3SJed Brown } 1073b1a0a8a3SJed Brown 1074b1a0a8a3SJed Brown #undef __FUNCT__ 1075f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMGetSubKSP" 1076b1a0a8a3SJed Brown /*@C 1077f746d493SDmitry Karpeev PCGASMGetSubKSP - Gets the local KSP contexts for all blocks on 1078b1a0a8a3SJed Brown this processor. 1079b1a0a8a3SJed Brown 1080b1a0a8a3SJed Brown Collective on PC iff first_local is requested 1081b1a0a8a3SJed Brown 1082b1a0a8a3SJed Brown Input Parameter: 1083b1a0a8a3SJed Brown . pc - the preconditioner context 1084b1a0a8a3SJed Brown 1085b1a0a8a3SJed Brown Output Parameters: 1086b1a0a8a3SJed Brown + n_local - the number of blocks on this processor or PETSC_NULL 1087b1a0a8a3SJed Brown . first_local - the global number of the first block on this processor or PETSC_NULL, 1088b1a0a8a3SJed Brown all processors must request or all must pass PETSC_NULL 1089b1a0a8a3SJed Brown - ksp - the array of KSP contexts 1090b1a0a8a3SJed Brown 1091b1a0a8a3SJed Brown Note: 1092f746d493SDmitry Karpeev After PCGASMGetSubKSP() the array of KSPes is not to be freed 1093b1a0a8a3SJed Brown 1094b1a0a8a3SJed Brown Currently for some matrix implementations only 1 block per processor 1095b1a0a8a3SJed Brown is supported. 1096b1a0a8a3SJed Brown 1097f746d493SDmitry Karpeev You must call KSPSetUp() before calling PCGASMGetSubKSP(). 1098b1a0a8a3SJed Brown 1099b1a0a8a3SJed Brown Level: advanced 1100b1a0a8a3SJed Brown 1101f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, get, sub, KSP, context 1102b1a0a8a3SJed Brown 11036a4f0f73SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetSubdomains(), PCGASMSetOverlap(), 1104f746d493SDmitry Karpeev PCGASMCreateSubdomains2D(), 1105b1a0a8a3SJed Brown @*/ 11067087cfbeSBarry Smith PetscErrorCode PCGASMGetSubKSP(PC pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp[]) 1107b1a0a8a3SJed Brown { 1108b1a0a8a3SJed Brown PetscErrorCode ierr; 1109b1a0a8a3SJed Brown 1110b1a0a8a3SJed Brown PetscFunctionBegin; 1111b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1112f746d493SDmitry Karpeev ierr = PetscUseMethod(pc,"PCGASMGetSubKSP_C",(PC,PetscInt*,PetscInt*,KSP **),(pc,n_local,first_local,ksp));CHKERRQ(ierr); 1113b1a0a8a3SJed Brown PetscFunctionReturn(0); 1114b1a0a8a3SJed Brown } 1115b1a0a8a3SJed Brown 1116b1a0a8a3SJed Brown /* -------------------------------------------------------------------------------------*/ 1117b1a0a8a3SJed Brown /*MC 1118f746d493SDmitry Karpeev PCGASM - Use the (restricted) additive Schwarz method, each block is (approximately) solved with 1119b1a0a8a3SJed Brown its own KSP object. 1120b1a0a8a3SJed Brown 1121b1a0a8a3SJed Brown Options Database Keys: 112206b43e7eSDmitry Karpeev + -pc_gasm_total_block_count <n> - Sets total number of local subdomains (known as blocks) to be distributed among processors 112306b43e7eSDmitry Karpeev . -pc_gasm_view_subdomains - activates the printing of subdomain indices in PCView(), -ksp_view or -snes_view 112406b43e7eSDmitry Karpeev . -pc_gasm_print_subdomains - activates the printing of subdomain indices in PCSetUp() 112506b43e7eSDmitry Karpeev . -pc_gasm_overlap <ovl> - Sets overlap by which to (automatically) extend local subdomains 1126f746d493SDmitry Karpeev - -pc_gasm_type [basic,restrict,interpolate,none] - Sets GASM type 1127b1a0a8a3SJed Brown 1128b1a0a8a3SJed Brown IMPORTANT: If you run with, for example, 3 blocks on 1 processor or 3 blocks on 3 processors you 1129f746d493SDmitry Karpeev will get a different convergence rate due to the default option of -pc_gasm_type restrict. Use 1130f746d493SDmitry Karpeev -pc_gasm_type basic to use the standard GASM. 1131b1a0a8a3SJed Brown 1132b1a0a8a3SJed Brown Notes: Each processor can have one or more blocks, but a block cannot be shared by more 1133b1a0a8a3SJed Brown than one processor. Defaults to one block per processor. 1134b1a0a8a3SJed Brown 1135b1a0a8a3SJed Brown To set options on the solvers for each block append -sub_ to all the KSP, and PC 1136b1a0a8a3SJed Brown options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly 1137b1a0a8a3SJed Brown 1138f746d493SDmitry Karpeev To set the options on the solvers separate for each block call PCGASMGetSubKSP() 1139b1a0a8a3SJed Brown and set the options directly on the resulting KSP object (you can access its PC 1140b1a0a8a3SJed Brown with KSPGetPC()) 1141b1a0a8a3SJed Brown 1142b1a0a8a3SJed Brown 1143b1a0a8a3SJed Brown Level: beginner 1144b1a0a8a3SJed Brown 1145b1a0a8a3SJed Brown Concepts: additive Schwarz method 1146b1a0a8a3SJed Brown 1147b1a0a8a3SJed Brown References: 1148b1a0a8a3SJed Brown An additive variant of the Schwarz alternating method for the case of many subregions 1149b1a0a8a3SJed Brown M Dryja, OB Widlund - Courant Institute, New York University Technical report 1150b1a0a8a3SJed Brown 1151b1a0a8a3SJed Brown Domain Decompositions: Parallel Multilevel Methods for Elliptic Partial Differential Equations, 1152b1a0a8a3SJed Brown Barry Smith, Petter Bjorstad, and William Gropp, Cambridge University Press, ISBN 0-521-49589-X. 1153b1a0a8a3SJed Brown 1154b1a0a8a3SJed Brown .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 115506b43e7eSDmitry Karpeev PCBJACOBI, PCGASMSetUseTrueLocal(), PCGASMGetSubKSP(), PCGASMSetSubdomains(), 11566a4f0f73SDmitry Karpeev PCGASMSetTotalSubdomains(), PCSetModifySubmatrices(), PCGASMSetOverlap(), PCGASMSetType() 1157b1a0a8a3SJed Brown 1158b1a0a8a3SJed Brown M*/ 1159b1a0a8a3SJed Brown 1160b1a0a8a3SJed Brown EXTERN_C_BEGIN 1161b1a0a8a3SJed Brown #undef __FUNCT__ 1162f746d493SDmitry Karpeev #define __FUNCT__ "PCCreate_GASM" 11637087cfbeSBarry Smith PetscErrorCode PCCreate_GASM(PC pc) 1164b1a0a8a3SJed Brown { 1165b1a0a8a3SJed Brown PetscErrorCode ierr; 1166f746d493SDmitry Karpeev PC_GASM *osm; 1167b1a0a8a3SJed Brown 1168b1a0a8a3SJed Brown PetscFunctionBegin; 1169f746d493SDmitry Karpeev ierr = PetscNewLog(pc,PC_GASM,&osm);CHKERRQ(ierr); 1170ab3e923fSDmitry Karpeev osm->N = PETSC_DECIDE; 117106b43e7eSDmitry Karpeev osm->n = PETSC_DECIDE; 1172ab3e923fSDmitry Karpeev osm->nmax = 0; 1173b1a0a8a3SJed Brown osm->overlap = 1; 1174b1a0a8a3SJed Brown osm->ksp = 0; 11756a4f0f73SDmitry Karpeev osm->gorestriction = 0; 11766a4f0f73SDmitry Karpeev osm->girestriction = 0; 1177ab3e923fSDmitry Karpeev osm->gx = 0; 1178ab3e923fSDmitry Karpeev osm->gy = 0; 1179b1a0a8a3SJed Brown osm->x = 0; 1180b1a0a8a3SJed Brown osm->y = 0; 11816a4f0f73SDmitry Karpeev osm->ois = 0; 11826a4f0f73SDmitry Karpeev osm->iis = 0; 1183b1a0a8a3SJed Brown osm->pmat = 0; 1184f746d493SDmitry Karpeev osm->type = PC_GASM_RESTRICT; 11856a4f0f73SDmitry Karpeev osm->same_subdomain_solvers = PETSC_TRUE; 1186b1a0a8a3SJed Brown osm->sort_indices = PETSC_TRUE; 1187b1a0a8a3SJed Brown 1188b1a0a8a3SJed Brown pc->data = (void*)osm; 1189f746d493SDmitry Karpeev pc->ops->apply = PCApply_GASM; 1190f746d493SDmitry Karpeev pc->ops->applytranspose = PCApplyTranspose_GASM; 1191f746d493SDmitry Karpeev pc->ops->setup = PCSetUp_GASM; 1192a06653b4SBarry Smith pc->ops->reset = PCReset_GASM; 1193f746d493SDmitry Karpeev pc->ops->destroy = PCDestroy_GASM; 1194f746d493SDmitry Karpeev pc->ops->setfromoptions = PCSetFromOptions_GASM; 1195f746d493SDmitry Karpeev pc->ops->setuponblocks = PCSetUpOnBlocks_GASM; 1196f746d493SDmitry Karpeev pc->ops->view = PCView_GASM; 1197b1a0a8a3SJed Brown pc->ops->applyrichardson = 0; 1198b1a0a8a3SJed Brown 119906b43e7eSDmitry Karpeev ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMSetSubdomains_C","PCGASMSetSubdomains_GASM", 120006b43e7eSDmitry Karpeev PCGASMSetSubdomains_GASM);CHKERRQ(ierr); 12016a4f0f73SDmitry Karpeev ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMSetTotalSubdomains_C","PCGASMSetTotalSubdomains_GASM", 12026a4f0f73SDmitry Karpeev PCGASMSetTotalSubdomains_GASM);CHKERRQ(ierr); 1203f746d493SDmitry Karpeev ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMSetOverlap_C","PCGASMSetOverlap_GASM", 1204f746d493SDmitry Karpeev PCGASMSetOverlap_GASM);CHKERRQ(ierr); 1205f746d493SDmitry Karpeev ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMSetType_C","PCGASMSetType_GASM", 1206f746d493SDmitry Karpeev PCGASMSetType_GASM);CHKERRQ(ierr); 1207f746d493SDmitry Karpeev ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMSetSortIndices_C","PCGASMSetSortIndices_GASM", 1208f746d493SDmitry Karpeev PCGASMSetSortIndices_GASM);CHKERRQ(ierr); 1209f746d493SDmitry Karpeev ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMGetSubKSP_C","PCGASMGetSubKSP_GASM", 1210f746d493SDmitry Karpeev PCGASMGetSubKSP_GASM);CHKERRQ(ierr); 1211b1a0a8a3SJed Brown PetscFunctionReturn(0); 1212b1a0a8a3SJed Brown } 1213b1a0a8a3SJed Brown EXTERN_C_END 1214b1a0a8a3SJed Brown 1215b1a0a8a3SJed Brown 1216b1a0a8a3SJed Brown #undef __FUNCT__ 121706b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMCreateLocalSubdomains" 1218b1a0a8a3SJed Brown /*@C 121906b43e7eSDmitry Karpeev PCGASMCreateLocalSubdomains - Creates n local index sets for the overlapping 122006b43e7eSDmitry Karpeev Schwarz preconditioner for a any problem based on its matrix. 1221b1a0a8a3SJed Brown 1222b1a0a8a3SJed Brown Collective 1223b1a0a8a3SJed Brown 1224b1a0a8a3SJed Brown Input Parameters: 1225b1a0a8a3SJed Brown + A - The global matrix operator 12266a4f0f73SDmitry Karpeev . overlap - amount of overlap in outer subdomains 122706b43e7eSDmitry Karpeev - n - the number of local subdomains 1228b1a0a8a3SJed Brown 1229b1a0a8a3SJed Brown Output Parameters: 12306a4f0f73SDmitry Karpeev + iis - the array of index sets defining the local inner subdomains (on which the correction is applied) 12316a4f0f73SDmitry Karpeev - ois - the array of index sets defining the local outer subdomains (on which the residual is computed) 1232b1a0a8a3SJed Brown 1233b1a0a8a3SJed Brown Level: advanced 1234b1a0a8a3SJed Brown 12356a4f0f73SDmitry Karpeev Note: this generates n nonoverlapping local inner subdomains on PETSC_COMM_SELF; 12366a4f0f73SDmitry Karpeev PCGASM will generate the overlap from these if you use them in PCGASMSetSubdomains() and set a 12376a4f0f73SDmitry Karpeev nonzero overlap with PCGASMSetOverlap() 1238b1a0a8a3SJed Brown 1239b1a0a8a3SJed Brown In the Fortran version you must provide the array outis[] already allocated of length n. 1240b1a0a8a3SJed Brown 1241f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, create, subdomains, unstructured grid 1242b1a0a8a3SJed Brown 124306b43e7eSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMDestroySubdomains() 1244b1a0a8a3SJed Brown @*/ 12456a4f0f73SDmitry Karpeev PetscErrorCode PCGASMCreateLocalSubdomains(Mat A, PetscInt overlap, PetscInt n, IS* iis[], IS* ois[]) 1246b1a0a8a3SJed Brown { 1247b1a0a8a3SJed Brown MatPartitioning mpart; 1248b1a0a8a3SJed Brown const char *prefix; 124911bd1e4dSLisandro Dalcin PetscErrorCode (*f)(Mat,MatReuse,Mat*); 1250b1a0a8a3SJed Brown PetscMPIInt size; 1251b1a0a8a3SJed Brown PetscInt i,j,rstart,rend,bs; 125211bd1e4dSLisandro Dalcin PetscBool isbaij = PETSC_FALSE,foundpart = PETSC_FALSE; 1253b1a0a8a3SJed Brown Mat Ad = PETSC_NULL, adj; 1254b1a0a8a3SJed Brown IS ispart,isnumb,*is; 1255b1a0a8a3SJed Brown PetscErrorCode ierr; 1256b1a0a8a3SJed Brown 1257b1a0a8a3SJed Brown PetscFunctionBegin; 1258b1a0a8a3SJed Brown PetscValidHeaderSpecific(A,MAT_CLASSID,1); 12596a4f0f73SDmitry Karpeev PetscValidPointer(iis,4); 1260b1a0a8a3SJed Brown if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of local blocks must be > 0, n = %D",n); 1261b1a0a8a3SJed Brown 1262b1a0a8a3SJed Brown /* Get prefix, row distribution, and block size */ 1263b1a0a8a3SJed Brown ierr = MatGetOptionsPrefix(A,&prefix);CHKERRQ(ierr); 1264b1a0a8a3SJed Brown ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 1265b1a0a8a3SJed Brown ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr); 1266b1a0a8a3SJed 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); 1267b1a0a8a3SJed Brown 1268b1a0a8a3SJed Brown /* Get diagonal block from matrix if possible */ 1269b1a0a8a3SJed Brown ierr = MPI_Comm_size(((PetscObject)A)->comm,&size);CHKERRQ(ierr); 1270b1a0a8a3SJed Brown ierr = PetscObjectQueryFunction((PetscObject)A,"MatGetDiagonalBlock_C",(void (**)(void))&f);CHKERRQ(ierr); 1271b1a0a8a3SJed Brown if (f) { 127211bd1e4dSLisandro Dalcin ierr = MatGetDiagonalBlock(A,&Ad);CHKERRQ(ierr); 1273b1a0a8a3SJed Brown } else if (size == 1) { 127411bd1e4dSLisandro Dalcin Ad = A; 1275b1a0a8a3SJed Brown } 1276b1a0a8a3SJed Brown if (Ad) { 1277251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)Ad,MATSEQBAIJ,&isbaij);CHKERRQ(ierr); 1278251f4c67SDmitry Karpeev if (!isbaij) {ierr = PetscObjectTypeCompare((PetscObject)Ad,MATSEQSBAIJ,&isbaij);CHKERRQ(ierr);} 1279b1a0a8a3SJed Brown } 1280b1a0a8a3SJed Brown if (Ad && n > 1) { 1281b1a0a8a3SJed Brown PetscBool match,done; 1282b1a0a8a3SJed Brown /* Try to setup a good matrix partitioning if available */ 1283b1a0a8a3SJed Brown ierr = MatPartitioningCreate(PETSC_COMM_SELF,&mpart);CHKERRQ(ierr); 1284b1a0a8a3SJed Brown ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr); 1285b1a0a8a3SJed Brown ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr); 1286251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGCURRENT,&match);CHKERRQ(ierr); 1287b1a0a8a3SJed Brown if (!match) { 1288251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGSQUARE,&match);CHKERRQ(ierr); 1289b1a0a8a3SJed Brown } 1290b1a0a8a3SJed Brown if (!match) { /* assume a "good" partitioner is available */ 1291b1a0a8a3SJed Brown PetscInt na,*ia,*ja; 1292b1a0a8a3SJed Brown ierr = MatGetRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr); 1293b1a0a8a3SJed Brown if (done) { 1294b1a0a8a3SJed Brown /* Build adjacency matrix by hand. Unfortunately a call to 1295b1a0a8a3SJed Brown MatConvert(Ad,MATMPIADJ,MAT_INITIAL_MATRIX,&adj) will 1296b1a0a8a3SJed Brown remove the block-aij structure and we cannot expect 1297b1a0a8a3SJed Brown MatPartitioning to split vertices as we need */ 1298b1a0a8a3SJed Brown PetscInt i,j,*row,len,nnz,cnt,*iia=0,*jja=0; 1299b1a0a8a3SJed Brown nnz = 0; 1300b1a0a8a3SJed Brown for (i=0; i<na; i++) { /* count number of nonzeros */ 1301b1a0a8a3SJed Brown len = ia[i+1] - ia[i]; 1302b1a0a8a3SJed Brown row = ja + ia[i]; 1303b1a0a8a3SJed Brown for (j=0; j<len; j++) { 1304b1a0a8a3SJed Brown if (row[j] == i) { /* don't count diagonal */ 1305b1a0a8a3SJed Brown len--; break; 1306b1a0a8a3SJed Brown } 1307b1a0a8a3SJed Brown } 1308b1a0a8a3SJed Brown nnz += len; 1309b1a0a8a3SJed Brown } 1310b1a0a8a3SJed Brown ierr = PetscMalloc((na+1)*sizeof(PetscInt),&iia);CHKERRQ(ierr); 1311b1a0a8a3SJed Brown ierr = PetscMalloc((nnz)*sizeof(PetscInt),&jja);CHKERRQ(ierr); 1312b1a0a8a3SJed Brown nnz = 0; 1313b1a0a8a3SJed Brown iia[0] = 0; 1314b1a0a8a3SJed Brown for (i=0; i<na; i++) { /* fill adjacency */ 1315b1a0a8a3SJed Brown cnt = 0; 1316b1a0a8a3SJed Brown len = ia[i+1] - ia[i]; 1317b1a0a8a3SJed Brown row = ja + ia[i]; 1318b1a0a8a3SJed Brown for (j=0; j<len; j++) { 1319b1a0a8a3SJed Brown if (row[j] != i) { /* if not diagonal */ 1320b1a0a8a3SJed Brown jja[nnz+cnt++] = row[j]; 1321b1a0a8a3SJed Brown } 1322b1a0a8a3SJed Brown } 1323b1a0a8a3SJed Brown nnz += cnt; 1324b1a0a8a3SJed Brown iia[i+1] = nnz; 1325b1a0a8a3SJed Brown } 1326b1a0a8a3SJed Brown /* Partitioning of the adjacency matrix */ 1327b1a0a8a3SJed Brown ierr = MatCreateMPIAdj(PETSC_COMM_SELF,na,na,iia,jja,PETSC_NULL,&adj);CHKERRQ(ierr); 1328b1a0a8a3SJed Brown ierr = MatPartitioningSetAdjacency(mpart,adj);CHKERRQ(ierr); 1329b1a0a8a3SJed Brown ierr = MatPartitioningSetNParts(mpart,n);CHKERRQ(ierr); 1330b1a0a8a3SJed Brown ierr = MatPartitioningApply(mpart,&ispart);CHKERRQ(ierr); 1331b1a0a8a3SJed Brown ierr = ISPartitioningToNumbering(ispart,&isnumb);CHKERRQ(ierr); 13326bf464f9SBarry Smith ierr = MatDestroy(&adj);CHKERRQ(ierr); 1333b1a0a8a3SJed Brown foundpart = PETSC_TRUE; 1334b1a0a8a3SJed Brown } 1335b1a0a8a3SJed Brown ierr = MatRestoreRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr); 1336b1a0a8a3SJed Brown } 13376bf464f9SBarry Smith ierr = MatPartitioningDestroy(&mpart);CHKERRQ(ierr); 1338b1a0a8a3SJed Brown } 1339b1a0a8a3SJed Brown ierr = PetscMalloc(n*sizeof(IS),&is);CHKERRQ(ierr); 1340b1a0a8a3SJed Brown if (!foundpart) { 1341b1a0a8a3SJed Brown 1342b1a0a8a3SJed Brown /* Partitioning by contiguous chunks of rows */ 1343b1a0a8a3SJed Brown 1344b1a0a8a3SJed Brown PetscInt mbs = (rend-rstart)/bs; 1345b1a0a8a3SJed Brown PetscInt start = rstart; 1346b1a0a8a3SJed Brown for (i=0; i<n; i++) { 1347b1a0a8a3SJed Brown PetscInt count = (mbs/n + ((mbs % n) > i)) * bs; 1348b1a0a8a3SJed Brown ierr = ISCreateStride(PETSC_COMM_SELF,count,start,1,&is[i]);CHKERRQ(ierr); 1349b1a0a8a3SJed Brown start += count; 1350b1a0a8a3SJed Brown } 1351b1a0a8a3SJed Brown 1352b1a0a8a3SJed Brown } else { 1353b1a0a8a3SJed Brown 1354b1a0a8a3SJed Brown /* Partitioning by adjacency of diagonal block */ 1355b1a0a8a3SJed Brown 1356b1a0a8a3SJed Brown const PetscInt *numbering; 1357b1a0a8a3SJed Brown PetscInt *count,nidx,*indices,*newidx,start=0; 1358b1a0a8a3SJed Brown /* Get node count in each partition */ 1359b1a0a8a3SJed Brown ierr = PetscMalloc(n*sizeof(PetscInt),&count);CHKERRQ(ierr); 1360b1a0a8a3SJed Brown ierr = ISPartitioningCount(ispart,n,count);CHKERRQ(ierr); 1361b1a0a8a3SJed Brown if (isbaij && bs > 1) { /* adjust for the block-aij case */ 1362b1a0a8a3SJed Brown for (i=0; i<n; i++) count[i] *= bs; 1363b1a0a8a3SJed Brown } 1364b1a0a8a3SJed Brown /* Build indices from node numbering */ 1365b1a0a8a3SJed Brown ierr = ISGetLocalSize(isnumb,&nidx);CHKERRQ(ierr); 1366b1a0a8a3SJed Brown ierr = PetscMalloc(nidx*sizeof(PetscInt),&indices);CHKERRQ(ierr); 1367b1a0a8a3SJed Brown for (i=0; i<nidx; i++) indices[i] = i; /* needs to be initialized */ 1368b1a0a8a3SJed Brown ierr = ISGetIndices(isnumb,&numbering);CHKERRQ(ierr); 1369b1a0a8a3SJed Brown ierr = PetscSortIntWithPermutation(nidx,numbering,indices);CHKERRQ(ierr); 1370b1a0a8a3SJed Brown ierr = ISRestoreIndices(isnumb,&numbering);CHKERRQ(ierr); 1371b1a0a8a3SJed Brown if (isbaij && bs > 1) { /* adjust for the block-aij case */ 1372b1a0a8a3SJed Brown ierr = PetscMalloc(nidx*bs*sizeof(PetscInt),&newidx);CHKERRQ(ierr); 1373b1a0a8a3SJed Brown for (i=0; i<nidx; i++) 1374b1a0a8a3SJed Brown for (j=0; j<bs; j++) 1375b1a0a8a3SJed Brown newidx[i*bs+j] = indices[i]*bs + j; 1376b1a0a8a3SJed Brown ierr = PetscFree(indices);CHKERRQ(ierr); 1377b1a0a8a3SJed Brown nidx *= bs; 1378b1a0a8a3SJed Brown indices = newidx; 1379b1a0a8a3SJed Brown } 1380b1a0a8a3SJed Brown /* Shift to get global indices */ 1381b1a0a8a3SJed Brown for (i=0; i<nidx; i++) indices[i] += rstart; 1382b1a0a8a3SJed Brown 1383b1a0a8a3SJed Brown /* Build the index sets for each block */ 1384b1a0a8a3SJed Brown for (i=0; i<n; i++) { 1385b1a0a8a3SJed Brown ierr = ISCreateGeneral(PETSC_COMM_SELF,count[i],&indices[start],PETSC_COPY_VALUES,&is[i]);CHKERRQ(ierr); 1386b1a0a8a3SJed Brown ierr = ISSort(is[i]);CHKERRQ(ierr); 1387b1a0a8a3SJed Brown start += count[i]; 1388b1a0a8a3SJed Brown } 1389b1a0a8a3SJed Brown 1390b1a0a8a3SJed Brown ierr = PetscFree(count); 1391b1a0a8a3SJed Brown ierr = PetscFree(indices); 1392fcfd50ebSBarry Smith ierr = ISDestroy(&isnumb);CHKERRQ(ierr); 1393fcfd50ebSBarry Smith ierr = ISDestroy(&ispart);CHKERRQ(ierr); 1394b1a0a8a3SJed Brown } 13956a4f0f73SDmitry Karpeev *iis = is; 13966a4f0f73SDmitry Karpeev if(!ois) PetscFunctionReturn(0); 13976a4f0f73SDmitry Karpeev /* 13986a4f0f73SDmitry Karpeev Initially make outer subdomains the same as inner subdomains. If nonzero additional overlap 13996a4f0f73SDmitry Karpeev has been requested, copy the inner subdomains over so they can be modified. 14006a4f0f73SDmitry Karpeev */ 14016a4f0f73SDmitry Karpeev ierr = PetscMalloc(n*sizeof(IS),ois);CHKERRQ(ierr); 14026a4f0f73SDmitry Karpeev for (i=0; i<n; ++i) { 14036a4f0f73SDmitry Karpeev if (overlap > 0) { /* With positive overlap, (*iis)[i] will be modified */ 14046a4f0f73SDmitry Karpeev ierr = ISDuplicate((*iis)[i],(*ois)+i);CHKERRQ(ierr); 14056a4f0f73SDmitry Karpeev ierr = ISCopy((*iis)[i],(*ois)[i]);CHKERRQ(ierr); 14066a4f0f73SDmitry Karpeev } else { 14076a4f0f73SDmitry Karpeev ierr = PetscObjectReference((PetscObject)(*iis)[i]);CHKERRQ(ierr); 14086a4f0f73SDmitry Karpeev (*ois)[i] = (*iis)[i]; 14096a4f0f73SDmitry Karpeev } 14106a4f0f73SDmitry Karpeev } 14116a4f0f73SDmitry Karpeev if (overlap > 0) { 14126a4f0f73SDmitry Karpeev /* Extend the "overlapping" regions by a number of steps */ 14136a4f0f73SDmitry Karpeev ierr = MatIncreaseOverlap(A,n,*ois,overlap);CHKERRQ(ierr); 14146a4f0f73SDmitry Karpeev } 1415b1a0a8a3SJed Brown PetscFunctionReturn(0); 1416b1a0a8a3SJed Brown } 1417b1a0a8a3SJed Brown 1418b1a0a8a3SJed Brown #undef __FUNCT__ 1419f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMDestroySubdomains" 1420b1a0a8a3SJed Brown /*@C 1421f746d493SDmitry Karpeev PCGASMDestroySubdomains - Destroys the index sets created with 142206b43e7eSDmitry Karpeev PCGASMCreateLocalSubdomains() or PCGASMCreateSubdomains2D. Should be 142306b43e7eSDmitry Karpeev called after setting subdomains with PCGASMSetSubdomains(). 1424b1a0a8a3SJed Brown 1425b1a0a8a3SJed Brown Collective 1426b1a0a8a3SJed Brown 1427b1a0a8a3SJed Brown Input Parameters: 1428b1a0a8a3SJed Brown + n - the number of index sets 14296a4f0f73SDmitry Karpeev . iis - the array of inner subdomains, 14306a4f0f73SDmitry Karpeev - ois - the array of outer subdomains, can be PETSC_NULL 1431b1a0a8a3SJed Brown 14326a4f0f73SDmitry Karpeev Level: intermediate 14336a4f0f73SDmitry Karpeev 14346a4f0f73SDmitry Karpeev Notes: this is merely a convenience subroutine that walks each list, 14356a4f0f73SDmitry Karpeev destroys each IS on the list, and then frees the list. 1436b1a0a8a3SJed Brown 1437f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, create, subdomains, unstructured grid 1438b1a0a8a3SJed Brown 143906b43e7eSDmitry Karpeev .seealso: PCGASMCreateLocalSubdomains(), PCGASMSetSubdomains() 1440b1a0a8a3SJed Brown @*/ 14416a4f0f73SDmitry Karpeev PetscErrorCode PCGASMDestroySubdomains(PetscInt n, IS iis[], IS ois[]) 1442b1a0a8a3SJed Brown { 1443b1a0a8a3SJed Brown PetscInt i; 1444b1a0a8a3SJed Brown PetscErrorCode ierr; 1445b1a0a8a3SJed Brown PetscFunctionBegin; 1446b1a0a8a3SJed Brown if (n <= 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"n must be > 0: n = %D",n); 14476a4f0f73SDmitry Karpeev PetscValidPointer(iis,2); 14486a4f0f73SDmitry Karpeev for (i=0; i<n; i++) { ierr = ISDestroy(&iis[i]);CHKERRQ(ierr); } 14496a4f0f73SDmitry Karpeev ierr = PetscFree(iis);CHKERRQ(ierr); 14506a4f0f73SDmitry Karpeev if (ois) { 14516a4f0f73SDmitry Karpeev for (i=0; i<n; i++) { ierr = ISDestroy(&ois[i]);CHKERRQ(ierr); } 14526a4f0f73SDmitry Karpeev ierr = PetscFree(ois);CHKERRQ(ierr); 1453b1a0a8a3SJed Brown } 1454b1a0a8a3SJed Brown PetscFunctionReturn(0); 1455b1a0a8a3SJed Brown } 1456b1a0a8a3SJed Brown 1457af538404SDmitry Karpeev 1458af538404SDmitry Karpeev #define PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,xleft_loc,ylow_loc,xright_loc,yhigh_loc,n) \ 1459af538404SDmitry Karpeev { \ 1460af538404SDmitry Karpeev PetscInt first_row = first/M, last_row = last/M+1; \ 1461af538404SDmitry Karpeev /* \ 1462af538404SDmitry Karpeev Compute ylow_loc and yhigh_loc so that (ylow_loc,xleft) and (yhigh_loc,xright) are the corners \ 1463af538404SDmitry Karpeev of the bounding box of the intersection of the subdomain with the local ownership range (local \ 1464af538404SDmitry Karpeev subdomain). \ 1465af538404SDmitry Karpeev Also compute xleft_loc and xright_loc as the lower and upper bounds on the first and last rows \ 1466af538404SDmitry Karpeev of the intersection. \ 1467af538404SDmitry Karpeev */ \ 1468af538404SDmitry Karpeev /* ylow_loc is the grid row containing the first element of the local sumbdomain */ \ 1469eec7e3faSDmitry Karpeev *ylow_loc = PetscMax(first_row,ylow); \ 1470af538404SDmitry Karpeev /* xleft_loc is the offset of first element of the local subdomain within its grid row (might actually be outside the local subdomain) */ \ 1471eec7e3faSDmitry Karpeev *xleft_loc = *ylow_loc==first_row?PetscMax(first%M,xleft):xleft; \ 1472af538404SDmitry Karpeev /* yhigh_loc is the grid row above the last local subdomain element */ \ 1473eec7e3faSDmitry Karpeev *yhigh_loc = PetscMin(last_row,yhigh); \ 1474af538404SDmitry Karpeev /* xright is the offset of the end of the local subdomain within its grid row (might actually be outside the local subdomain) */ \ 1475eec7e3faSDmitry Karpeev *xright_loc = *yhigh_loc==last_row?PetscMin(xright,last%M):xright; \ 1476af538404SDmitry Karpeev /* Now compute the size of the local subdomain n. */ \ 1477c3518bceSDmitry Karpeev *n = 0; \ 1478eec7e3faSDmitry Karpeev if(*ylow_loc < *yhigh_loc) { \ 1479af538404SDmitry Karpeev PetscInt width = xright-xleft; \ 1480eec7e3faSDmitry Karpeev *n += width*(*yhigh_loc-*ylow_loc-1); \ 1481eec7e3faSDmitry Karpeev *n += PetscMin(PetscMax(*xright_loc-xleft,0),width); \ 1482eec7e3faSDmitry Karpeev *n -= PetscMin(PetscMax(*xleft_loc-xleft,0), width); \ 1483af538404SDmitry Karpeev }\ 1484af538404SDmitry Karpeev } 1485af538404SDmitry Karpeev 1486c3518bceSDmitry Karpeev 1487eec7e3faSDmitry Karpeev 1488eec7e3faSDmitry Karpeev #undef __FUNCT__ 1489f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMCreateSubdomains2D" 1490b1a0a8a3SJed Brown /*@ 1491f746d493SDmitry Karpeev PCGASMCreateSubdomains2D - Creates the index sets for the overlapping Schwarz 1492b1a0a8a3SJed Brown preconditioner for a two-dimensional problem on a regular grid. 1493b1a0a8a3SJed Brown 1494af538404SDmitry Karpeev Collective 1495b1a0a8a3SJed Brown 1496b1a0a8a3SJed Brown Input Parameters: 149706b43e7eSDmitry Karpeev + M, N - the global number of grid points in the x and y directions 1498af538404SDmitry Karpeev . Mdomains, Ndomains - the global number of subdomains in the x and y directions 1499b1a0a8a3SJed Brown . dof - degrees of freedom per node 1500b1a0a8a3SJed Brown - overlap - overlap in mesh lines 1501b1a0a8a3SJed Brown 1502b1a0a8a3SJed Brown Output Parameters: 1503af538404SDmitry Karpeev + Nsub - the number of local subdomains created 15046a4f0f73SDmitry Karpeev . iis - array of index sets defining inner (nonoverlapping) subdomains 15056a4f0f73SDmitry Karpeev - ois - array of index sets defining outer (overlapping, if overlap > 0) subdomains 1506b1a0a8a3SJed Brown 1507b1a0a8a3SJed Brown 1508b1a0a8a3SJed Brown Level: advanced 1509b1a0a8a3SJed Brown 1510f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, create, subdomains, 2D, regular grid 1511b1a0a8a3SJed Brown 15126a4f0f73SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetSubdomains(), PCGASMGetSubKSP(), 1513f746d493SDmitry Karpeev PCGASMSetOverlap() 1514b1a0a8a3SJed Brown @*/ 15156a4f0f73SDmitry Karpeev PetscErrorCode PCGASMCreateSubdomains2D(PC pc, PetscInt M,PetscInt N,PetscInt Mdomains,PetscInt Ndomains,PetscInt dof,PetscInt overlap, PetscInt *nsub,IS **iis,IS **ois) 1516b1a0a8a3SJed Brown { 1517b1a0a8a3SJed Brown PetscErrorCode ierr; 15182bbb417fSDmitry Karpeev PetscMPIInt size, rank; 15192bbb417fSDmitry Karpeev PetscInt i, j; 15202bbb417fSDmitry Karpeev PetscInt maxheight, maxwidth; 15212bbb417fSDmitry Karpeev PetscInt xstart, xleft, xright, xleft_loc, xright_loc; 15222bbb417fSDmitry Karpeev PetscInt ystart, ylow, yhigh, ylow_loc, yhigh_loc; 1523eec7e3faSDmitry Karpeev PetscInt x[2][2], y[2][2], n[2]; 15242bbb417fSDmitry Karpeev PetscInt first, last; 15252bbb417fSDmitry Karpeev PetscInt nidx, *idx; 15262bbb417fSDmitry Karpeev PetscInt ii,jj,s,q,d; 1527af538404SDmitry Karpeev PetscInt k,kk; 15282bbb417fSDmitry Karpeev PetscMPIInt color; 15292bbb417fSDmitry Karpeev MPI_Comm comm, subcomm; 15306a4f0f73SDmitry Karpeev IS **xis = 0, **is = ois, **is_local = iis; 1531d34fcf5fSBarry Smith 1532b1a0a8a3SJed Brown PetscFunctionBegin; 15332bbb417fSDmitry Karpeev ierr = PetscObjectGetComm((PetscObject)pc, &comm);CHKERRQ(ierr); 15342bbb417fSDmitry Karpeev ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 15352bbb417fSDmitry Karpeev ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 15362bbb417fSDmitry Karpeev ierr = MatGetOwnershipRange(pc->pmat, &first, &last);CHKERRQ(ierr); 1537d34fcf5fSBarry 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) " 15382bbb417fSDmitry Karpeev "does not respect the number of degrees of freedom per grid point %D", first, last, dof); 1539d34fcf5fSBarry Smith 1540af538404SDmitry Karpeev /* Determine the number of domains with nonzero intersections with the local ownership range. */ 15412bbb417fSDmitry Karpeev s = 0; 1542af538404SDmitry Karpeev ystart = 0; 1543af538404SDmitry Karpeev for (j=0; j<Ndomains; ++j) { 1544af538404SDmitry Karpeev maxheight = N/Ndomains + ((N % Ndomains) > j); /* Maximal height of subdomain */ 1545af538404SDmitry 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); 1546eec7e3faSDmitry Karpeev /* Vertical domain limits with an overlap. */ 1547eec7e3faSDmitry Karpeev ylow = PetscMax(ystart - overlap,0); 1548eec7e3faSDmitry Karpeev yhigh = PetscMin(ystart + maxheight + overlap,N); 1549b1a0a8a3SJed Brown xstart = 0; 1550af538404SDmitry Karpeev for (i=0; i<Mdomains; ++i) { 1551af538404SDmitry Karpeev maxwidth = M/Mdomains + ((M % Mdomains) > i); /* Maximal width of subdomain */ 1552af538404SDmitry 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); 1553eec7e3faSDmitry Karpeev /* Horizontal domain limits with an overlap. */ 1554eec7e3faSDmitry Karpeev xleft = PetscMax(xstart - overlap,0); 1555eec7e3faSDmitry Karpeev xright = PetscMin(xstart + maxwidth + overlap,M); 1556af538404SDmitry Karpeev /* 1557af538404SDmitry Karpeev Determine whether this subdomain intersects this processor's ownership range of pc->pmat. 1558af538404SDmitry Karpeev */ 1559c3518bceSDmitry Karpeev PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,(&xleft_loc),(&ylow_loc),(&xright_loc),(&yhigh_loc),(&nidx)); 1560af538404SDmitry Karpeev if(nidx) { 1561af538404SDmitry Karpeev ++s; 1562af538404SDmitry Karpeev } 1563af538404SDmitry Karpeev xstart += maxwidth; 1564af538404SDmitry Karpeev }/* for(i = 0; i < Mdomains; ++i) */ 1565af538404SDmitry Karpeev ystart += maxheight; 1566af538404SDmitry Karpeev }/* for(j = 0; j < Ndomains; ++j) */ 1567af538404SDmitry Karpeev /* Now we can allocate the necessary number of ISs. */ 1568af538404SDmitry Karpeev *nsub = s; 1569af538404SDmitry Karpeev ierr = PetscMalloc((*nsub)*sizeof(IS*),is);CHKERRQ(ierr); 1570af538404SDmitry Karpeev ierr = PetscMalloc((*nsub)*sizeof(IS*),is_local);CHKERRQ(ierr); 1571af538404SDmitry Karpeev s = 0; 1572af538404SDmitry Karpeev ystart = 0; 1573af538404SDmitry Karpeev for (j=0; j<Ndomains; ++j) { 1574af538404SDmitry Karpeev maxheight = N/Ndomains + ((N % Ndomains) > j); /* Maximal height of subdomain */ 1575af538404SDmitry 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); 1576af538404SDmitry Karpeev /* Vertical domain limits with an overlap. */ 1577eec7e3faSDmitry Karpeev y[0][0] = PetscMax(ystart - overlap,0); 1578eec7e3faSDmitry Karpeev y[0][1] = PetscMin(ystart + maxheight + overlap,N); 1579eec7e3faSDmitry Karpeev /* Vertical domain limits without an overlap. */ 1580eec7e3faSDmitry Karpeev y[1][0] = ystart; 1581eec7e3faSDmitry Karpeev y[1][1] = PetscMin(ystart + maxheight,N); 1582eec7e3faSDmitry Karpeev xstart = 0; 1583af538404SDmitry Karpeev for (i=0; i<Mdomains; ++i) { 1584af538404SDmitry Karpeev maxwidth = M/Mdomains + ((M % Mdomains) > i); /* Maximal width of subdomain */ 1585af538404SDmitry 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); 1586af538404SDmitry Karpeev /* Horizontal domain limits with an overlap. */ 1587eec7e3faSDmitry Karpeev x[0][0] = PetscMax(xstart - overlap,0); 1588eec7e3faSDmitry Karpeev x[0][1] = PetscMin(xstart + maxwidth + overlap,M); 1589eec7e3faSDmitry Karpeev /* Horizontal domain limits without an overlap. */ 1590eec7e3faSDmitry Karpeev x[1][0] = xstart; 1591eec7e3faSDmitry Karpeev x[1][1] = PetscMin(xstart+maxwidth,M); 15922bbb417fSDmitry Karpeev /* 15932bbb417fSDmitry Karpeev Determine whether this domain intersects this processor's ownership range of pc->pmat. 15942bbb417fSDmitry Karpeev Do this twice: first for the domains with overlaps, and once without. 15952bbb417fSDmitry Karpeev During the first pass create the subcommunicators, and use them on the second pass as well. 15962bbb417fSDmitry Karpeev */ 15972bbb417fSDmitry Karpeev for(q = 0; q < 2; ++q) { 15982bbb417fSDmitry Karpeev /* 15992bbb417fSDmitry Karpeev domain limits, (xleft, xright) and (ylow, yheigh) are adjusted 16002bbb417fSDmitry Karpeev according to whether the domain with an overlap or without is considered. 16012bbb417fSDmitry Karpeev */ 16022bbb417fSDmitry Karpeev xleft = x[q][0]; xright = x[q][1]; 16032bbb417fSDmitry Karpeev ylow = y[q][0]; yhigh = y[q][1]; 1604c3518bceSDmitry Karpeev PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,(&xleft_loc),(&ylow_loc),(&xright_loc),(&yhigh_loc),(&nidx)); 1605af538404SDmitry Karpeev nidx *= dof; 1606eec7e3faSDmitry Karpeev n[q] = nidx; 16072bbb417fSDmitry Karpeev /* 1608eec7e3faSDmitry Karpeev Based on the counted number of indices in the local domain *with an overlap*, 1609af538404SDmitry Karpeev construct a subcommunicator of all the processors supporting this domain. 1610eec7e3faSDmitry Karpeev Observe that a domain with an overlap might have nontrivial local support, 1611eec7e3faSDmitry Karpeev while the domain without an overlap might not. Hence, the decision to participate 1612eec7e3faSDmitry Karpeev in the subcommunicator must be based on the domain with an overlap. 16132bbb417fSDmitry Karpeev */ 16142bbb417fSDmitry Karpeev if (q == 0) { 16152bbb417fSDmitry Karpeev if(nidx) { 16162bbb417fSDmitry Karpeev color = 1; 1617d34fcf5fSBarry Smith } else { 16182bbb417fSDmitry Karpeev color = MPI_UNDEFINED; 16192bbb417fSDmitry Karpeev } 16202bbb417fSDmitry Karpeev ierr = MPI_Comm_split(comm, color, rank, &subcomm);CHKERRQ(ierr); 16212bbb417fSDmitry Karpeev } 1622af538404SDmitry Karpeev /* 1623eec7e3faSDmitry Karpeev Proceed only if the number of local indices *with an overlap* is nonzero. 1624af538404SDmitry Karpeev */ 1625eec7e3faSDmitry Karpeev if (n[0]) { 1626af538404SDmitry Karpeev if(q == 0) { 16276a4f0f73SDmitry Karpeev xis = is; 1628af538404SDmitry Karpeev } 1629af538404SDmitry Karpeev if (q == 1) { 1630af538404SDmitry Karpeev /* 1631eec7e3faSDmitry Karpeev The IS for the no-overlap subdomain shares a communicator with the overlapping domain. 1632af538404SDmitry Karpeev Moreover, if the overlap is zero, the two ISs are identical. 1633af538404SDmitry Karpeev */ 1634b1a0a8a3SJed Brown if (overlap == 0) { 1635eec7e3faSDmitry Karpeev (*is_local)[s] = (*is)[s]; 16362bbb417fSDmitry Karpeev ierr = PetscObjectReference((PetscObject)(*is)[s]);CHKERRQ(ierr); 16372bbb417fSDmitry Karpeev continue; 1638d34fcf5fSBarry Smith } else { 16396a4f0f73SDmitry Karpeev xis = is_local; 1640eec7e3faSDmitry Karpeev subcomm = ((PetscObject)(*is)[s])->comm; 16412bbb417fSDmitry Karpeev } 1642af538404SDmitry Karpeev }/* if(q == 1) */ 1643eec7e3faSDmitry Karpeev idx = PETSC_NULL; 16442bbb417fSDmitry Karpeev ierr = PetscMalloc(nidx*sizeof(PetscInt),&idx);CHKERRQ(ierr); 1645eec7e3faSDmitry Karpeev if(nidx) { 16462bbb417fSDmitry Karpeev k = 0; 16472bbb417fSDmitry Karpeev for (jj=ylow_loc; jj<yhigh_loc; ++jj) { 1648af538404SDmitry Karpeev PetscInt x0 = (jj==ylow_loc)?xleft_loc:xleft; 1649af538404SDmitry Karpeev PetscInt x1 = (jj==yhigh_loc-1)?xright_loc:xright; 1650af538404SDmitry Karpeev kk = dof*(M*jj + x0); 16512bbb417fSDmitry Karpeev for (ii=x0; ii<x1; ++ii) { 16522bbb417fSDmitry Karpeev for(d = 0; d < dof; ++d) { 16532bbb417fSDmitry Karpeev idx[k++] = kk++; 1654b1a0a8a3SJed Brown } 1655b1a0a8a3SJed Brown } 1656b1a0a8a3SJed Brown } 1657eec7e3faSDmitry Karpeev } 16586a4f0f73SDmitry Karpeev ierr = ISCreateGeneral(subcomm,nidx,idx,PETSC_OWN_POINTER,(*xis)+s);CHKERRQ(ierr); 1659eec7e3faSDmitry Karpeev }/* if(n[0]) */ 16602bbb417fSDmitry Karpeev }/* for(q = 0; q < 2; ++q) */ 1661eec7e3faSDmitry Karpeev if(n[0]) { 16622bbb417fSDmitry Karpeev ++s; 1663b1a0a8a3SJed Brown } 1664af538404SDmitry Karpeev xstart += maxwidth; 1665af538404SDmitry Karpeev }/* for(i = 0; i < Mdomains; ++i) */ 16662bbb417fSDmitry Karpeev ystart += maxheight; 1667af538404SDmitry Karpeev }/* for(j = 0; j < Ndomains; ++j) */ 1668b1a0a8a3SJed Brown PetscFunctionReturn(0); 1669b1a0a8a3SJed Brown } 1670b1a0a8a3SJed Brown 1671b1a0a8a3SJed Brown #undef __FUNCT__ 167206b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMGetSubdomains" 1673b1a0a8a3SJed Brown /*@C 167406b43e7eSDmitry Karpeev PCGASMGetSubdomains - Gets the subdomains supported on this processor 167506b43e7eSDmitry Karpeev for the additive Schwarz preconditioner. 1676b1a0a8a3SJed Brown 1677b1a0a8a3SJed Brown Not Collective 1678b1a0a8a3SJed Brown 1679b1a0a8a3SJed Brown Input Parameter: 1680b1a0a8a3SJed Brown . pc - the preconditioner context 1681b1a0a8a3SJed Brown 1682b1a0a8a3SJed Brown Output Parameters: 1683b1a0a8a3SJed Brown + n - the number of subdomains for this processor (default value = 1) 16846a4f0f73SDmitry Karpeev . iis - the index sets that define the inner subdomains (without overlap) supported on this processor (can be PETSC_NULL) 16856a4f0f73SDmitry Karpeev - ois - the index sets that define the outer subdomains (with overlap) supported on this processor (can be PETSC_NULL) 1686b1a0a8a3SJed Brown 1687b1a0a8a3SJed Brown 1688b1a0a8a3SJed Brown Notes: 16896a4f0f73SDmitry Karpeev The user is responsible for destroying the ISs and freeing the returned arrays. 1690b1a0a8a3SJed Brown The IS numbering is in the parallel, global numbering of the vector. 1691b1a0a8a3SJed Brown 1692b1a0a8a3SJed Brown Level: advanced 1693b1a0a8a3SJed Brown 169406b43e7eSDmitry Karpeev .keywords: PC, GASM, get, subdomains, additive Schwarz 1695b1a0a8a3SJed Brown 16966a4f0f73SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetOverlap(), PCGASMGetSubKSP(), 169706b43e7eSDmitry Karpeev PCGASMCreateSubdomains2D(), PCGASMSetSubdomains(), PCGASMGetSubmatrices() 1698b1a0a8a3SJed Brown @*/ 16996a4f0f73SDmitry Karpeev PetscErrorCode PCGASMGetSubdomains(PC pc,PetscInt *n,IS *iis[],IS *ois[]) 1700b1a0a8a3SJed Brown { 1701f746d493SDmitry Karpeev PC_GASM *osm; 1702b1a0a8a3SJed Brown PetscErrorCode ierr; 1703b1a0a8a3SJed Brown PetscBool match; 17046a4f0f73SDmitry Karpeev PetscInt i; 1705b1a0a8a3SJed Brown PetscFunctionBegin; 1706b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1707251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match);CHKERRQ(ierr); 17086a4f0f73SDmitry Karpeev if (!match) 17096a4f0f73SDmitry Karpeev SETERRQ2(((PetscObject)pc)->comm, PETSC_ERR_ARG_WRONG, "Incorrect object type: expected %s, got %s instead", PCGASM, ((PetscObject)pc)->type_name); 1710f746d493SDmitry Karpeev osm = (PC_GASM*)pc->data; 1711ab3e923fSDmitry Karpeev if (n) *n = osm->n; 17126a4f0f73SDmitry Karpeev if(iis) { 17136a4f0f73SDmitry Karpeev ierr = PetscMalloc(osm->n*sizeof(IS), iis); CHKERRQ(ierr); 17146a4f0f73SDmitry Karpeev } 17156a4f0f73SDmitry Karpeev if(ois) { 17166a4f0f73SDmitry Karpeev ierr = PetscMalloc(osm->n*sizeof(IS), ois); CHKERRQ(ierr); 17176a4f0f73SDmitry Karpeev } 17186a4f0f73SDmitry Karpeev if(iis || ois) { 17196a4f0f73SDmitry Karpeev for(i = 0; i < osm->n; ++i) { 17206a4f0f73SDmitry Karpeev if(iis) (*iis)[i] = osm->iis[i]; 17216a4f0f73SDmitry Karpeev if(ois) (*ois)[i] = osm->ois[i]; 17226a4f0f73SDmitry Karpeev } 1723b1a0a8a3SJed Brown } 1724b1a0a8a3SJed Brown PetscFunctionReturn(0); 1725b1a0a8a3SJed Brown } 1726b1a0a8a3SJed Brown 1727b1a0a8a3SJed Brown #undef __FUNCT__ 172806b43e7eSDmitry Karpeev #define __FUNCT__ "PCGASMGetSubmatrices" 1729b1a0a8a3SJed Brown /*@C 173006b43e7eSDmitry Karpeev PCGASMGetSubmatrices - Gets the local submatrices (for this processor 1731b1a0a8a3SJed Brown only) for the additive Schwarz preconditioner. 1732b1a0a8a3SJed Brown 1733b1a0a8a3SJed Brown Not Collective 1734b1a0a8a3SJed Brown 1735b1a0a8a3SJed Brown Input Parameter: 1736b1a0a8a3SJed Brown . pc - the preconditioner context 1737b1a0a8a3SJed Brown 1738b1a0a8a3SJed Brown Output Parameters: 1739b1a0a8a3SJed Brown + n - the number of matrices for this processor (default value = 1) 1740b1a0a8a3SJed Brown - mat - the matrices 1741b1a0a8a3SJed Brown 174206b43e7eSDmitry Karpeev Notes: matrices returned by this routine have the same communicators as the index sets (IS) 174306b43e7eSDmitry Karpeev used to define subdomains in PCGASMSetSubdomains(), or PETSC_COMM_SELF, if the 17446a4f0f73SDmitry Karpeev subdomains were defined using PCGASMSetTotalSubdomains(). 1745b1a0a8a3SJed Brown Level: advanced 1746b1a0a8a3SJed Brown 1747f746d493SDmitry Karpeev .keywords: PC, GASM, set, local, subdomains, additive Schwarz, block Jacobi 1748b1a0a8a3SJed Brown 17496a4f0f73SDmitry Karpeev .seealso: PCGASMSetTotalSubdomain(), PCGASMSetOverlap(), PCGASMGetSubKSP(), 175006b43e7eSDmitry Karpeev PCGASMCreateSubdomains2D(), PCGASMSetSubdomains(), PCGASMGetSubdomains() 1751b1a0a8a3SJed Brown @*/ 175206b43e7eSDmitry Karpeev PetscErrorCode PCGASMGetSubmatrices(PC pc,PetscInt *n,Mat *mat[]) 1753b1a0a8a3SJed Brown { 1754f746d493SDmitry Karpeev PC_GASM *osm; 1755b1a0a8a3SJed Brown PetscErrorCode ierr; 1756b1a0a8a3SJed Brown PetscBool match; 1757b1a0a8a3SJed Brown 1758b1a0a8a3SJed Brown PetscFunctionBegin; 1759b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1760b1a0a8a3SJed Brown PetscValidIntPointer(n,2); 1761b1a0a8a3SJed Brown if (mat) PetscValidPointer(mat,3); 1762b1a0a8a3SJed Brown if (!pc->setupcalled) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must call after KSPSetUP() or PCSetUp()."); 1763251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match);CHKERRQ(ierr); 176406b43e7eSDmitry Karpeev if (!match) SETERRQ2(((PetscObject)pc)->comm, PETSC_ERR_ARG_WRONG, "Expected %s, got %s instead", PCGASM, ((PetscObject)pc)->type_name); 1765f746d493SDmitry Karpeev osm = (PC_GASM*)pc->data; 1766ab3e923fSDmitry Karpeev if (n) *n = osm->n; 1767b1a0a8a3SJed Brown if (mat) *mat = osm->pmat; 176806b43e7eSDmitry Karpeev 1769b1a0a8a3SJed Brown PetscFunctionReturn(0); 1770b1a0a8a3SJed Brown } 1771