1b1a0a8a3SJed Brown 2b1a0a8a3SJed Brown /* 3f746d493SDmitry Karpeev This file defines an "generalized" additive Schwarz preconditioner for any Mat implementation. 444fe18e5SDmitry Karpeev In this version each processor may have any number of subdomains and a subdomain may live on multiple 5f746d493SDmitry Karpeev processors. 6b1a0a8a3SJed Brown 7af538404SDmitry Karpeev N - total number of subdomains on all processors 8f746d493SDmitry Karpeev n - actual number of subdomains on this processor 9f746d493SDmitry Karpeev nmax - maximum number of subdomains per processor 10b1a0a8a3SJed Brown */ 11c6db04a5SJed Brown #include <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. */ 19f746d493SDmitry Karpeev IS gis, gis_local; /* merged ISs */ 20f746d493SDmitry Karpeev VecScatter grestriction; /* merged restriction */ 21f746d493SDmitry Karpeev VecScatter gprolongation; /* merged prolongation */ 22b1a0a8a3SJed Brown IS *is; /* index set that defines each overlapping subdomain */ 23f746d493SDmitry Karpeev IS *is_local; /* index set that defines each local subdomain (same as subdomain with the overlap removed); may be NULL */ 24f746d493SDmitry Karpeev Mat *pmat; /* subdomain block matrices */ 25f746d493SDmitry Karpeev PCGASMType type; /* use reduced interpolation, restriction or both */ 26b1a0a8a3SJed Brown PetscBool type_set; /* if user set this value (so won't change it for symmetric problems) */ 27b1a0a8a3SJed Brown PetscBool same_local_solves; /* 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__ 32af538404SDmitry Karpeev #define __FUNCT__ "PCGASMPrintSubdomains" 33af538404SDmitry Karpeev static PetscErrorCode PCGASMPrintSubdomains(PC pc) 34af538404SDmitry Karpeev { 35af538404SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 36af538404SDmitry Karpeev const char *prefix; 37af538404SDmitry Karpeev char fname[PETSC_MAX_PATH_LEN+1]; 38af538404SDmitry Karpeev PetscViewer viewer, sviewer; 39af538404SDmitry Karpeev PetscInt i,j,nidx; 40af538404SDmitry Karpeev const PetscInt *idx; 41af538404SDmitry Karpeev char *cidx; 42af538404SDmitry Karpeev PetscBool found; 43eec7e3faSDmitry Karpeev PetscMPIInt size, rank; 44af538404SDmitry Karpeev PetscErrorCode ierr; 45af538404SDmitry Karpeev 46af538404SDmitry Karpeev PetscFunctionBegin; 47eec7e3faSDmitry Karpeev ierr = MPI_Comm_size(((PetscObject)pc)->comm, &size);CHKERRQ(ierr); 48eec7e3faSDmitry Karpeev ierr = MPI_Comm_rank(((PetscObject)pc)->comm, &rank);CHKERRQ(ierr); 49af538404SDmitry Karpeev ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 50af538404SDmitry Karpeev ierr = PetscOptionsGetString(prefix,"-pc_gasm_print_subdomains",fname,PETSC_MAX_PATH_LEN,&found);CHKERRQ(ierr); 51af538404SDmitry Karpeev if (!found) { ierr = PetscStrcpy(fname,"stdout");CHKERRQ(ierr); }; 52af538404SDmitry Karpeev for (i=0;i<osm->n;++i) { 53af538404SDmitry Karpeev ierr = PetscViewerASCIIOpen(((PetscObject)((osm->is)[i]))->comm,fname,&viewer);CHKERRQ(ierr); 54af538404SDmitry Karpeev ierr = ISGetLocalSize(osm->is[i], &nidx);CHKERRQ(ierr); 55eec7e3faSDmitry Karpeev /* 56eec7e3faSDmitry Karpeev No more than 15 characters per index plus a space. 57eec7e3faSDmitry Karpeev PetscViewerStringSPrintf requires a string of size at least 2, so use (nidx+1) instead of nidx, 58eec7e3faSDmitry Karpeev in case nidx == 0. That will take care of the space for the trailing '\0' as well. 59eec7e3faSDmitry Karpeev For nidx == 0, the whole string 16 '\0'. 60eec7e3faSDmitry Karpeev */ 61eec7e3faSDmitry Karpeev ierr = PetscMalloc(sizeof(char)*(16*(nidx+1)+1), &cidx);CHKERRQ(ierr); 62af538404SDmitry Karpeev ierr = ISGetIndices(osm->is[i], &idx);CHKERRQ(ierr); 63eec7e3faSDmitry Karpeev ierr = PetscViewerStringOpen(((PetscObject)(osm->is[i]))->comm, cidx, 16*(nidx+1), &sviewer);CHKERRQ(ierr); 64af538404SDmitry Karpeev for(j = 0; j < nidx; ++j) { 65af538404SDmitry Karpeev ierr = PetscViewerStringSPrintf(sviewer, "%D ", idx[j]);CHKERRQ(ierr); 66af538404SDmitry Karpeev } 67af538404SDmitry Karpeev ierr = PetscViewerDestroy(sviewer);CHKERRQ(ierr); 68af538404SDmitry Karpeev ierr = ISRestoreIndices(osm->is[i],&idx);CHKERRQ(ierr); 69eec7e3faSDmitry Karpeev 70*7b23a99aSBarry Smith ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr); 71af538404SDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer, "Subdomain with overlap\n");CHKERRQ(ierr); 72af538404SDmitry Karpeev ierr = PetscViewerASCIISynchronizedPrintf(viewer, "%s", cidx);CHKERRQ(ierr); 73af538404SDmitry Karpeev ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 74*7b23a99aSBarry Smith ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr); 75af538404SDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer, "\n");CHKERRQ(ierr); 76af538404SDmitry Karpeev ierr = PetscFree(cidx);CHKERRQ(ierr); 77af538404SDmitry Karpeev ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr); 78af538404SDmitry Karpeev if (osm->is_local) { 79af538404SDmitry Karpeev ierr = PetscViewerASCIIOpen(((PetscObject)((osm->is)[i]))->comm,fname,&viewer);CHKERRQ(ierr); 80af538404SDmitry Karpeev ierr = ISGetLocalSize(osm->is_local[i], &nidx);CHKERRQ(ierr); 81eec7e3faSDmitry Karpeev /* 82eec7e3faSDmitry Karpeev No more than 15 characters per index plus a space. 83eec7e3faSDmitry Karpeev PetscViewerStringSPrintf requires a string of size at least 2, so use (nidx+1) instead of nidx, 84eec7e3faSDmitry Karpeev in case nidx == 0. That will take care of the space for the trailing '\0' as well. 85eec7e3faSDmitry Karpeev For nidx == 0, the whole string 16 '\0'. 86eec7e3faSDmitry Karpeev */ 87eec7e3faSDmitry Karpeev ierr = PetscMalloc(sizeof(char)*(16*(nidx+1)+1), &cidx);CHKERRQ(ierr); 88af538404SDmitry Karpeev ierr = ISGetIndices(osm->is_local[i], &idx);CHKERRQ(ierr); 89eec7e3faSDmitry Karpeev ierr = PetscViewerStringOpen(((PetscObject)(osm->is_local[i]))->comm, cidx, 16*(nidx+1), &sviewer);CHKERRQ(ierr); 90af538404SDmitry Karpeev for(j = 0; j < nidx; ++j) { 91af538404SDmitry Karpeev ierr = PetscViewerStringSPrintf(sviewer, "%D ", idx[j]);CHKERRQ(ierr); 92af538404SDmitry Karpeev } 93af538404SDmitry Karpeev ierr = PetscViewerDestroy(sviewer);CHKERRQ(ierr); 94af538404SDmitry Karpeev ierr = ISRestoreIndices(osm->is_local[i],&idx);CHKERRQ(ierr); 95eec7e3faSDmitry Karpeev 96*7b23a99aSBarry Smith ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr); 97af538404SDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer, "Subdomain without overlap\n");CHKERRQ(ierr); 98af538404SDmitry Karpeev ierr = PetscViewerASCIISynchronizedPrintf(viewer, "%s", cidx);CHKERRQ(ierr); 99af538404SDmitry Karpeev ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 100*7b23a99aSBarry Smith ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr); 101af538404SDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer, "\n");CHKERRQ(ierr); 102af538404SDmitry Karpeev ierr = PetscFree(cidx);CHKERRQ(ierr); 103af538404SDmitry Karpeev ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr); 104af538404SDmitry Karpeev } 105af538404SDmitry Karpeev } 106af538404SDmitry Karpeev PetscFunctionReturn(0); 107af538404SDmitry Karpeev } 108af538404SDmitry Karpeev 109af538404SDmitry Karpeev 110af538404SDmitry Karpeev #undef __FUNCT__ 111f746d493SDmitry Karpeev #define __FUNCT__ "PCView_GASM" 112f746d493SDmitry Karpeev static PetscErrorCode PCView_GASM(PC pc,PetscViewer viewer) 113b1a0a8a3SJed Brown { 114f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 115af538404SDmitry Karpeev const char *prefix; 116b1a0a8a3SJed Brown PetscErrorCode ierr; 117af538404SDmitry Karpeev PetscMPIInt rank, size; 118b1a0a8a3SJed Brown PetscInt i,bsz; 119af538404SDmitry Karpeev PetscBool iascii,isstring, print_subdomains; 120b1a0a8a3SJed Brown PetscViewer sviewer; 121b1a0a8a3SJed Brown 122b1a0a8a3SJed Brown 123b1a0a8a3SJed Brown PetscFunctionBegin; 124af538404SDmitry Karpeev ierr = MPI_Comm_size(((PetscObject)pc)->comm, &size);CHKERRQ(ierr); 125af538404SDmitry Karpeev ierr = MPI_Comm_rank(((PetscObject)pc)->comm, &rank);CHKERRQ(ierr); 126b1a0a8a3SJed Brown ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 127b1a0a8a3SJed Brown ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr); 128af538404SDmitry Karpeev ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 129af538404SDmitry Karpeev ierr = PetscOptionsGetBool(prefix,"-pc_gasm_print_subdomains",&print_subdomains,PETSC_NULL);CHKERRQ(ierr); 130b1a0a8a3SJed Brown if (iascii) { 131af538404SDmitry Karpeev char overlaps[256] = "user-defined overlap",subdomains[256] = "total subdomains set"; 132b1a0a8a3SJed Brown if (osm->overlap >= 0) {ierr = PetscSNPrintf(overlaps,sizeof overlaps,"amount of overlap = %D",osm->overlap);CHKERRQ(ierr);} 133af538404SDmitry Karpeev if (osm->nmax > 0) {ierr = PetscSNPrintf(subdomains,sizeof subdomains,"max number of local subdomains = %D",osm->nmax);CHKERRQ(ierr);} 134*7b23a99aSBarry Smith ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_TRUE);CHKERRQ(ierr); 135eec7e3faSDmitry Karpeev ierr = PetscViewerASCIISynchronizedPrintf(viewer," [%d:%d] number of locally-supported subdomains = %D\n",(int)rank,(int)size,osm->n);CHKERRQ(ierr); 136af538404SDmitry Karpeev ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 137*7b23a99aSBarry Smith ierr = PetscViewerASCIISynchronizedAllow(viewer,PETSC_FALSE);CHKERRQ(ierr); 138af538404SDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer," Generalized additive Schwarz: %s, %s\n",subdomains,overlaps);CHKERRQ(ierr); 139f746d493SDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer," Generalized additive Schwarz: restriction/interpolation type - %s\n",PCGASMTypes[osm->type]);CHKERRQ(ierr); 140b1a0a8a3SJed Brown ierr = MPI_Comm_rank(((PetscObject)pc)->comm,&rank);CHKERRQ(ierr); 141b1a0a8a3SJed Brown if (osm->same_local_solves) { 142b1a0a8a3SJed Brown if (osm->ksp) { 143af538404SDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer," Local solve is same for all subdomains, in the following KSP and PC objects:\n");CHKERRQ(ierr); 144b1a0a8a3SJed Brown ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr); 145b1a0a8a3SJed Brown if (!rank) { 146b1a0a8a3SJed Brown ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 147b1a0a8a3SJed Brown ierr = KSPView(osm->ksp[0],sviewer);CHKERRQ(ierr); 148b1a0a8a3SJed Brown ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 149b1a0a8a3SJed Brown } 150b1a0a8a3SJed Brown ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr); 151b1a0a8a3SJed Brown } 152b1a0a8a3SJed Brown } else { 153b1a0a8a3SJed Brown ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 154af538404SDmitry Karpeev ierr = PetscViewerASCIIPrintf(viewer," Local solve info for each subdomain is in the following KSP and PC objects:\n");CHKERRQ(ierr); 155b1a0a8a3SJed Brown ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 156b1a0a8a3SJed Brown ierr = PetscViewerASCIIPrintf(viewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr); 157f746d493SDmitry Karpeev for (i=0; i<osm->nmax; i++) { 158b1a0a8a3SJed Brown ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr); 159f746d493SDmitry Karpeev if (i < osm->n) { 160b1a0a8a3SJed Brown ierr = ISGetLocalSize(osm->is[i],&bsz);CHKERRQ(ierr); 161*7b23a99aSBarry Smith ierr = PetscViewerASCIISynchronizedAllow(sviewer,PETSC_TRUE);CHKERRQ(ierr); 162af538404SDmitry Karpeev ierr = PetscViewerASCIISynchronizedPrintf(sviewer,"[%d:%d] local subdomain number %D, size = %D\n",(int)rank,(int)size,i,bsz);CHKERRQ(ierr); 163b1a0a8a3SJed Brown ierr = KSPView(osm->ksp[i],sviewer);CHKERRQ(ierr); 164b1a0a8a3SJed Brown ierr = PetscViewerASCIISynchronizedPrintf(sviewer,"- - - - - - - - - - - - - - - - - -\n");CHKERRQ(ierr); 165*7b23a99aSBarry Smith ierr = PetscViewerFlush(sviewer);CHKERRQ(ierr); 166*7b23a99aSBarry Smith ierr = PetscViewerASCIISynchronizedAllow(sviewer,PETSC_FALSE);CHKERRQ(ierr); 167b1a0a8a3SJed Brown } 168b1a0a8a3SJed Brown ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr); 169b1a0a8a3SJed Brown } 170b1a0a8a3SJed Brown ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 171b1a0a8a3SJed Brown ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 172b1a0a8a3SJed Brown } 173b1a0a8a3SJed Brown } else if (isstring) { 174af538404SDmitry Karpeev ierr = PetscViewerStringSPrintf(viewer," subdomains=%D, overlap=%D, type=%s",osm->n,osm->overlap,PCGASMTypes[osm->type]);CHKERRQ(ierr); 175b1a0a8a3SJed Brown ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr); 176b1a0a8a3SJed Brown if (osm->ksp) {ierr = KSPView(osm->ksp[0],sviewer);CHKERRQ(ierr);} 177b1a0a8a3SJed Brown ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr); 178b1a0a8a3SJed Brown } else { 179f746d493SDmitry Karpeev SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PCGASM",((PetscObject)viewer)->type_name); 180b1a0a8a3SJed Brown } 181b1a0a8a3SJed Brown ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 182af538404SDmitry Karpeev if(print_subdomains) { 183af538404SDmitry Karpeev ierr = PCGASMPrintSubdomains(pc);CHKERRQ(ierr); 184af538404SDmitry Karpeev } 185b1a0a8a3SJed Brown PetscFunctionReturn(0); 186b1a0a8a3SJed Brown } 187b1a0a8a3SJed Brown 188af538404SDmitry Karpeev 189af538404SDmitry Karpeev 190af538404SDmitry Karpeev 191af538404SDmitry Karpeev 192b1a0a8a3SJed Brown #undef __FUNCT__ 193f746d493SDmitry Karpeev #define __FUNCT__ "PCSetUp_GASM" 194f746d493SDmitry Karpeev static PetscErrorCode PCSetUp_GASM(PC pc) 195b1a0a8a3SJed Brown { 196f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 197b1a0a8a3SJed Brown PetscErrorCode ierr; 198b1a0a8a3SJed Brown PetscBool symset,flg; 19988389755SJed Brown PetscInt i; 200c10bc1a1SDmitry Karpeev PetscMPIInt rank, size; 201b1a0a8a3SJed Brown MatReuse scall = MAT_REUSE_MATRIX; 202b1a0a8a3SJed Brown KSP ksp; 203b1a0a8a3SJed Brown PC subpc; 204b1a0a8a3SJed Brown const char *prefix,*pprefix; 205f746d493SDmitry Karpeev PetscInt dn; /* Number of indices in a single subdomain assigned to this processor. */ 206f746d493SDmitry Karpeev const PetscInt *didx; /* Indices from a single subdomain assigned to this processor. */ 207f746d493SDmitry Karpeev PetscInt ddn; /* Number of indices in all subdomains assigned to this processor. */ 208f746d493SDmitry Karpeev PetscInt *ddidx; /* Indices of all subdomains assigned to this processor. */ 209f746d493SDmitry Karpeev IS gid; /* Identity IS of the size of all subdomains assigned to this processor. */ 210f746d493SDmitry Karpeev Vec x,y; 211f746d493SDmitry Karpeev PetscScalar *gxarray, *gyarray; 2126d98aee9SDmitry Karpeev PetscInt gfirst, glast; 213b1a0a8a3SJed Brown 214b1a0a8a3SJed Brown PetscFunctionBegin; 215c10bc1a1SDmitry Karpeev ierr = MPI_Comm_size(((PetscObject)pc)->comm,&size);CHKERRQ(ierr); 216c10bc1a1SDmitry Karpeev ierr = MPI_Comm_rank(((PetscObject)pc)->comm,&rank);CHKERRQ(ierr); 217b1a0a8a3SJed Brown if (!pc->setupcalled) { 218b1a0a8a3SJed Brown 219b1a0a8a3SJed Brown if (!osm->type_set) { 220b1a0a8a3SJed Brown ierr = MatIsSymmetricKnown(pc->pmat,&symset,&flg);CHKERRQ(ierr); 221f746d493SDmitry Karpeev if (symset && flg) { osm->type = PC_GASM_BASIC; } 222b1a0a8a3SJed Brown } 223b1a0a8a3SJed Brown 22444fe18e5SDmitry Karpeev if (osm->N == PETSC_DECIDE && osm->n < 1) { 225b1a0a8a3SJed Brown /* no subdomains given, use one per processor */ 22644fe18e5SDmitry Karpeev osm->nmax = osm->n = 1; 227b1a0a8a3SJed Brown ierr = MPI_Comm_size(((PetscObject)pc)->comm,&size);CHKERRQ(ierr); 228f746d493SDmitry Karpeev osm->N = size; 229f746d493SDmitry Karpeev } else if (osm->N == PETSC_DECIDE) { 230b1a0a8a3SJed Brown PetscInt inwork[2], outwork[2]; 231f746d493SDmitry Karpeev /* determine global number of subdomains and the max number of local subdomains */ 232f746d493SDmitry Karpeev inwork[0] = inwork[1] = osm->n; 233b1a0a8a3SJed Brown ierr = MPI_Allreduce(inwork,outwork,1,MPIU_2INT,PetscMaxSum_Op,((PetscObject)pc)->comm);CHKERRQ(ierr); 234f746d493SDmitry Karpeev osm->nmax = outwork[0]; 235f746d493SDmitry Karpeev osm->N = outwork[1]; 236b1a0a8a3SJed Brown } 237b1a0a8a3SJed Brown if (!osm->is){ /* create the index sets */ 238f746d493SDmitry Karpeev ierr = PCGASMCreateSubdomains(pc->pmat,osm->n,&osm->is);CHKERRQ(ierr); 239b1a0a8a3SJed Brown } 24093c1ec32SDmitry Karpeev if (!osm->is_local) { 24193c1ec32SDmitry Karpeev /* 24293c1ec32SDmitry Karpeev This indicates that osm->is should define a nonoverlapping decomposition 24393c1ec32SDmitry Karpeev (there is no way to really guarantee that if subdomains are set by the user through PCGASMSetLocalSubdomains, 24493c1ec32SDmitry Karpeev but the assumption is that either the user does the right thing, or subdomains in ossm->is have been created 24593c1ec32SDmitry Karpeev via PCGASMCreateSubdomains, which guarantees a nonoverlapping decomposition). 24693c1ec32SDmitry Karpeev Therefore, osm->is will be used to define osm->is_local. 24793c1ec32SDmitry Karpeev If a nonzero overlap has been requested by the user, then osm->is will be expanded and will overlap, 24893c1ec32SDmitry Karpeev so osm->is_local should obtain a copy of osm->is while they are still (presumably) nonoverlapping. 24993c1ec32SDmitry Karpeev Otherwise (no overlap has been requested), osm->is_local are simply aliases for osm->is. 25093c1ec32SDmitry Karpeev */ 251f746d493SDmitry Karpeev ierr = PetscMalloc(osm->n*sizeof(IS),&osm->is_local);CHKERRQ(ierr); 252f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 253b1a0a8a3SJed Brown if (osm->overlap > 0) { /* With positive overlap, osm->is[i] will be modified */ 254b1a0a8a3SJed Brown ierr = ISDuplicate(osm->is[i],&osm->is_local[i]);CHKERRQ(ierr); 255b1a0a8a3SJed Brown ierr = ISCopy(osm->is[i],osm->is_local[i]);CHKERRQ(ierr); 256b1a0a8a3SJed Brown } else { 257b1a0a8a3SJed Brown ierr = PetscObjectReference((PetscObject)osm->is[i]);CHKERRQ(ierr); 258b1a0a8a3SJed Brown osm->is_local[i] = osm->is[i]; 259b1a0a8a3SJed Brown } 260b1a0a8a3SJed Brown } 261b1a0a8a3SJed Brown } 26293c1ec32SDmitry Karpeev /* Beyond this point osm->is_local is not null. */ 26393c1ec32SDmitry Karpeev if (osm->overlap > 0) { 26493c1ec32SDmitry Karpeev /* Extend the "overlapping" regions by a number of steps */ 26593c1ec32SDmitry Karpeev ierr = MatIncreaseOverlap(pc->pmat,osm->n,osm->is,osm->overlap);CHKERRQ(ierr); 26693c1ec32SDmitry Karpeev } 267b1a0a8a3SJed Brown ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 268b1a0a8a3SJed Brown flg = PETSC_FALSE; 269f746d493SDmitry Karpeev ierr = PetscOptionsGetBool(prefix,"-pc_gasm_print_subdomains",&flg,PETSC_NULL);CHKERRQ(ierr); 270f746d493SDmitry Karpeev if (flg) { ierr = PCGASMPrintSubdomains(pc);CHKERRQ(ierr); } 271b1a0a8a3SJed Brown 272b1a0a8a3SJed Brown if (osm->sort_indices) { 273f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 274b1a0a8a3SJed Brown ierr = ISSort(osm->is[i]);CHKERRQ(ierr); 275b1a0a8a3SJed Brown ierr = ISSort(osm->is_local[i]);CHKERRQ(ierr); 276b1a0a8a3SJed Brown } 277b1a0a8a3SJed Brown } 278f746d493SDmitry Karpeev /* Merge the ISs, create merged vectors and scatter contexts. */ 27993c1ec32SDmitry Karpeev /* Restriction ISs. */ 280f746d493SDmitry Karpeev ddn = 0; 281f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 282f746d493SDmitry Karpeev ierr = ISGetLocalSize(osm->is[i],&dn);CHKERRQ(ierr); 283f746d493SDmitry Karpeev ddn += dn; 284f746d493SDmitry Karpeev } 285f746d493SDmitry Karpeev ierr = PetscMalloc(ddn*sizeof(PetscInt), &ddidx);CHKERRQ(ierr); 286f746d493SDmitry Karpeev ddn = 0; 287f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 288f746d493SDmitry Karpeev ierr = ISGetLocalSize(osm->is[i],&dn);CHKERRQ(ierr); 289f746d493SDmitry Karpeev ierr = ISGetIndices(osm->is[i],&didx);CHKERRQ(ierr); 290f746d493SDmitry Karpeev ierr = PetscMemcpy(ddidx+ddn, didx, sizeof(PetscInt)*dn);CHKERRQ(ierr); 291f746d493SDmitry Karpeev ierr = ISRestoreIndices(osm->is[i], &didx);CHKERRQ(ierr); 29293c1ec32SDmitry Karpeev ddn += dn; 293f746d493SDmitry Karpeev } 294f746d493SDmitry Karpeev ierr = ISCreateGeneral(((PetscObject)(pc))->comm, ddn, ddidx, PETSC_OWN_POINTER, &osm->gis);CHKERRQ(ierr); 295f746d493SDmitry Karpeev ierr = MatGetVecs(pc->pmat,&x,&y);CHKERRQ(ierr); 296f746d493SDmitry Karpeev ierr = VecCreateMPI(((PetscObject)pc)->comm, ddn, PETSC_DECIDE, &osm->gx);CHKERRQ(ierr); 297f746d493SDmitry Karpeev ierr = VecDuplicate(osm->gx,&osm->gy);CHKERRQ(ierr); 2986d98aee9SDmitry Karpeev ierr = VecGetOwnershipRange(osm->gx, &gfirst, &glast);CHKERRQ(ierr); 2996d98aee9SDmitry Karpeev ierr = ISCreateStride(((PetscObject)pc)->comm,ddn,gfirst,1, &gid);CHKERRQ(ierr); 300f746d493SDmitry Karpeev ierr = VecScatterCreate(x,osm->gis,osm->gx,gid, &(osm->grestriction));CHKERRQ(ierr); 301f746d493SDmitry Karpeev ierr = ISDestroy(gid);CHKERRQ(ierr); 30293c1ec32SDmitry Karpeev /* Prolongation ISs */ 30381a5c4aaSDmitry Karpeev { PetscInt dn_local; /* Number of indices in the local part of a single domain assigned to this processor. */ 304f746d493SDmitry Karpeev const PetscInt *didx_local; /* Global indices from the local part of a single domain assigned to this processor. */ 305f746d493SDmitry Karpeev PetscInt ddn_local; /* Number of indices in the local part of the disjoint union all domains assigned to this processor. */ 306f746d493SDmitry Karpeev PetscInt *ddidx_local; /* Global indices of the local part of the disjoint union of all domains assigned to this processor. */ 307f746d493SDmitry Karpeev /**/ 30893c1ec32SDmitry Karpeev ISLocalToGlobalMapping ltog; /* Map from global to local indices on the disjoint union of subdomains: "local" ind's run from 0 to ddn-1. */ 30993c1ec32SDmitry Karpeev PetscInt *ddidx_llocal; /* Mapped local indices of the disjoint union of local parts of subdomains. */ 310f746d493SDmitry Karpeev PetscInt ddn_llocal; /* Number of indices in ddidx_llocal; must equal ddn_local, or else gis_local is not a sub-IS of gis. */ 311f746d493SDmitry Karpeev IS gis_llocal; /* IS with ddidx_llocal indices. */ 312f746d493SDmitry Karpeev PetscInt j; 313f746d493SDmitry Karpeev ddn_local = 0; 314f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 315f746d493SDmitry Karpeev ierr = ISGetLocalSize(osm->is_local[i],&dn_local);CHKERRQ(ierr); 316f746d493SDmitry Karpeev ddn_local += dn_local; 317f746d493SDmitry Karpeev } 318f746d493SDmitry Karpeev ierr = PetscMalloc(ddn_local*sizeof(PetscInt), &ddidx_local);CHKERRQ(ierr); 319f746d493SDmitry Karpeev ddn_local = 0; 320f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 321f746d493SDmitry Karpeev ierr = ISGetLocalSize(osm->is_local[i],&dn_local);CHKERRQ(ierr); 322f746d493SDmitry Karpeev ierr = ISGetIndices(osm->is_local[i],&didx_local);CHKERRQ(ierr); 323f746d493SDmitry Karpeev ierr = PetscMemcpy(ddidx_local+ddn_local, didx_local, sizeof(PetscInt)*dn_local);CHKERRQ(ierr); 324f746d493SDmitry Karpeev ierr = ISRestoreIndices(osm->is_local[i], &didx_local);CHKERRQ(ierr); 325f746d493SDmitry Karpeev ddn_local += dn_local; 326f746d493SDmitry Karpeev } 327ab3e923fSDmitry Karpeev ierr = PetscMalloc(sizeof(PetscInt)*ddn_local, &ddidx_llocal);CHKERRQ(ierr); 328f746d493SDmitry Karpeev ierr = ISCreateGeneral(((PetscObject)pc)->comm, ddn_local, ddidx_local, PETSC_OWN_POINTER, &(osm->gis_local));CHKERRQ(ierr); 329f746d493SDmitry Karpeev ierr = ISLocalToGlobalMappingCreateIS(osm->gis,<og);CHKERRQ(ierr); 33081a5c4aaSDmitry Karpeev ierr = ISGlobalToLocalMappingApply(ltog,IS_GTOLM_DROP,ddn_local,ddidx_local,&ddn_llocal,ddidx_llocal);CHKERRQ(ierr); 331b1a0a8a3SJed Brown ierr = ISLocalToGlobalMappingDestroy(ltog);CHKERRQ(ierr); 33293c1ec32SDmitry Karpeev if (ddn_llocal != ddn_local) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_PLIB,"gis_local contains %D indices outside of gis", ddn_llocal - ddn_local); 333f746d493SDmitry Karpeev /* Now convert these localized indices into the global indices into the merged output vector. */ 3346d98aee9SDmitry Karpeev ierr = VecGetOwnershipRange(osm->gy, &gfirst, &glast);CHKERRQ(ierr); 335f746d493SDmitry Karpeev for(j=0; j < ddn_llocal; ++j) { 3366d98aee9SDmitry Karpeev ddidx_llocal[j] += gfirst; 337f746d493SDmitry Karpeev } 338f746d493SDmitry Karpeev ierr = ISCreateGeneral(PETSC_COMM_SELF,ddn_llocal,ddidx_llocal,PETSC_OWN_POINTER,&gis_llocal);CHKERRQ(ierr); 339f746d493SDmitry Karpeev ierr = VecScatterCreate(y,osm->gis_local,osm->gy,gis_llocal,&osm->gprolongation);CHKERRQ(ierr); 340f746d493SDmitry Karpeev ierr = ISDestroy(gis_llocal);CHKERRQ(ierr); 341b1a0a8a3SJed Brown } 342f746d493SDmitry Karpeev /* Create the subdomain work vectors. */ 343f746d493SDmitry Karpeev ierr = PetscMalloc(osm->n*sizeof(Vec),&osm->x);CHKERRQ(ierr); 344f746d493SDmitry Karpeev ierr = PetscMalloc(osm->n*sizeof(Vec),&osm->y);CHKERRQ(ierr); 3456d98aee9SDmitry Karpeev ierr = VecGetOwnershipRange(osm->gx, &gfirst, &glast);CHKERRQ(ierr); 346f746d493SDmitry Karpeev ierr = VecGetArray(osm->gx, &gxarray);CHKERRQ(ierr); 347f746d493SDmitry Karpeev ierr = VecGetArray(osm->gy, &gyarray);CHKERRQ(ierr); 348f746d493SDmitry Karpeev for (i=0, ddn=0; i<osm->n; ++i, ddn += dn) { 34986b96d36SDmitry Karpeev PetscInt dN; 350f746d493SDmitry Karpeev ierr = ISGetLocalSize(osm->is[i],&dn);CHKERRQ(ierr); 35186b96d36SDmitry Karpeev ierr = ISGetSize(osm->is[i],&dN);CHKERRQ(ierr); 35286b96d36SDmitry Karpeev ierr = VecCreateMPIWithArray(((PetscObject)(osm->is[i]))->comm,dn,dN,gxarray+ddn,&osm->x[i]);CHKERRQ(ierr); 35386b96d36SDmitry Karpeev ierr = VecCreateMPIWithArray(((PetscObject)(osm->is[i]))->comm,dn,dN,gyarray+ddn,&osm->y[i]);CHKERRQ(ierr); 354b1a0a8a3SJed Brown } 355f746d493SDmitry Karpeev ierr = VecRestoreArray(osm->gx, &gxarray);CHKERRQ(ierr); 356f746d493SDmitry Karpeev ierr = VecRestoreArray(osm->gy, &gyarray);CHKERRQ(ierr); 357f746d493SDmitry Karpeev ierr = VecDestroy(x);CHKERRQ(ierr); 358f746d493SDmitry Karpeev ierr = VecDestroy(y);CHKERRQ(ierr); 359b1a0a8a3SJed Brown /* Create the local solvers */ 360f746d493SDmitry Karpeev ierr = PetscMalloc(osm->n*sizeof(KSP *),&osm->ksp);CHKERRQ(ierr); 36144fe18e5SDmitry Karpeev for (i=0; i<osm->n; i++) { /* KSPs are local */ 36286b96d36SDmitry Karpeev ierr = KSPCreate(((PetscObject)(osm->is[i]))->comm,&ksp);CHKERRQ(ierr); 363b1a0a8a3SJed Brown ierr = PetscLogObjectParent(pc,ksp);CHKERRQ(ierr); 364b1a0a8a3SJed Brown ierr = PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,1);CHKERRQ(ierr); 365b1a0a8a3SJed Brown ierr = KSPSetType(ksp,KSPPREONLY);CHKERRQ(ierr); 366b1a0a8a3SJed Brown ierr = KSPGetPC(ksp,&subpc);CHKERRQ(ierr); 367b1a0a8a3SJed Brown ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 368b1a0a8a3SJed Brown ierr = KSPSetOptionsPrefix(ksp,prefix);CHKERRQ(ierr); 369b1a0a8a3SJed Brown ierr = KSPAppendOptionsPrefix(ksp,"sub_");CHKERRQ(ierr); 370b1a0a8a3SJed Brown osm->ksp[i] = ksp; 371b1a0a8a3SJed Brown } 372b1a0a8a3SJed Brown scall = MAT_INITIAL_MATRIX; 373b1a0a8a3SJed Brown 374b1a0a8a3SJed Brown } else { 375b1a0a8a3SJed Brown /* 376b1a0a8a3SJed Brown Destroy the blocks from the previous iteration 377b1a0a8a3SJed Brown */ 378b1a0a8a3SJed Brown if (pc->flag == DIFFERENT_NONZERO_PATTERN) { 379f746d493SDmitry Karpeev ierr = MatDestroyMatrices(osm->n,&osm->pmat);CHKERRQ(ierr); 380b1a0a8a3SJed Brown scall = MAT_INITIAL_MATRIX; 381b1a0a8a3SJed Brown } 382b1a0a8a3SJed Brown } 383b1a0a8a3SJed Brown 384b1a0a8a3SJed Brown /* 385f746d493SDmitry Karpeev Extract out the submatrices. 386b1a0a8a3SJed Brown */ 38781a5c4aaSDmitry Karpeev if(size > 1) { 3886d42056aSDmitry Karpeev ierr = MatGetSubMatricesParallel(pc->pmat,osm->n,osm->is, osm->is,scall,&osm->pmat);CHKERRQ(ierr); 38981a5c4aaSDmitry Karpeev } 39081a5c4aaSDmitry Karpeev else { 39181a5c4aaSDmitry Karpeev ierr = MatGetSubMatrices(pc->pmat,osm->n,osm->is, osm->is,scall,&osm->pmat);CHKERRQ(ierr); 39281a5c4aaSDmitry Karpeev } 393b1a0a8a3SJed Brown if (scall == MAT_INITIAL_MATRIX) { 394b1a0a8a3SJed Brown ierr = PetscObjectGetOptionsPrefix((PetscObject)pc->pmat,&pprefix);CHKERRQ(ierr); 395f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 396b1a0a8a3SJed Brown ierr = PetscLogObjectParent(pc,osm->pmat[i]);CHKERRQ(ierr); 397b1a0a8a3SJed Brown ierr = PetscObjectSetOptionsPrefix((PetscObject)osm->pmat[i],pprefix);CHKERRQ(ierr); 398b1a0a8a3SJed Brown } 399b1a0a8a3SJed Brown } 400b1a0a8a3SJed Brown 401b1a0a8a3SJed Brown /* Return control to the user so that the submatrices can be modified (e.g., to apply 402b1a0a8a3SJed Brown different boundary conditions for the submatrices than for the global problem) */ 403c10bc1a1SDmitry Karpeev ierr = PCModifySubMatrices(pc,osm->n,osm->is,osm->is,osm->pmat,pc->modifysubmatricesP);CHKERRQ(ierr); 404b1a0a8a3SJed Brown 405b1a0a8a3SJed Brown /* 406f746d493SDmitry Karpeev Loop over submatrices putting them into local ksp 407b1a0a8a3SJed Brown */ 408f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 409b1a0a8a3SJed Brown ierr = KSPSetOperators(osm->ksp[i],osm->pmat[i],osm->pmat[i],pc->flag);CHKERRQ(ierr); 410b1a0a8a3SJed Brown if (!pc->setupcalled) { 411b1a0a8a3SJed Brown ierr = KSPSetFromOptions(osm->ksp[i]);CHKERRQ(ierr); 412b1a0a8a3SJed Brown } 413b1a0a8a3SJed Brown } 414b1a0a8a3SJed Brown 415b1a0a8a3SJed Brown PetscFunctionReturn(0); 416b1a0a8a3SJed Brown } 417b1a0a8a3SJed Brown 418b1a0a8a3SJed Brown #undef __FUNCT__ 419f746d493SDmitry Karpeev #define __FUNCT__ "PCSetUpOnBlocks_GASM" 420f746d493SDmitry Karpeev static PetscErrorCode PCSetUpOnBlocks_GASM(PC pc) 421b1a0a8a3SJed Brown { 422f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 423b1a0a8a3SJed Brown PetscErrorCode ierr; 424b1a0a8a3SJed Brown PetscInt i; 425b1a0a8a3SJed Brown 426b1a0a8a3SJed Brown PetscFunctionBegin; 427f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 428b1a0a8a3SJed Brown ierr = KSPSetUp(osm->ksp[i]);CHKERRQ(ierr); 429b1a0a8a3SJed Brown } 430b1a0a8a3SJed Brown PetscFunctionReturn(0); 431b1a0a8a3SJed Brown } 432b1a0a8a3SJed Brown 433b1a0a8a3SJed Brown #undef __FUNCT__ 434f746d493SDmitry Karpeev #define __FUNCT__ "PCApply_GASM" 435f746d493SDmitry Karpeev static PetscErrorCode PCApply_GASM(PC pc,Vec x,Vec y) 436b1a0a8a3SJed Brown { 437f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 438b1a0a8a3SJed Brown PetscErrorCode ierr; 439f746d493SDmitry Karpeev PetscInt i; 440b1a0a8a3SJed Brown ScatterMode forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE; 441b1a0a8a3SJed Brown 442b1a0a8a3SJed Brown PetscFunctionBegin; 443b1a0a8a3SJed Brown /* 444b1a0a8a3SJed Brown Support for limiting the restriction or interpolation to only local 445b1a0a8a3SJed Brown subdomain values (leaving the other values 0). 446b1a0a8a3SJed Brown */ 447f746d493SDmitry Karpeev if (!(osm->type & PC_GASM_RESTRICT)) { 448b1a0a8a3SJed Brown forward = SCATTER_FORWARD_LOCAL; 449b1a0a8a3SJed Brown /* have to zero the work RHS since scatter may leave some slots empty */ 450f746d493SDmitry Karpeev ierr = VecZeroEntries(osm->gx);CHKERRQ(ierr); 451b1a0a8a3SJed Brown } 452f746d493SDmitry Karpeev if (!(osm->type & PC_GASM_INTERPOLATE)) { 453b1a0a8a3SJed Brown reverse = SCATTER_REVERSE_LOCAL; 454b1a0a8a3SJed Brown } 455b1a0a8a3SJed Brown 456f746d493SDmitry Karpeev ierr = VecScatterBegin(osm->grestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 457b1a0a8a3SJed Brown ierr = VecZeroEntries(y);CHKERRQ(ierr); 458f746d493SDmitry Karpeev ierr = VecScatterEnd(osm->grestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 45986b96d36SDmitry Karpeev /* do the subdomain solves */ 46086b96d36SDmitry Karpeev for (i=0; i<osm->n; ++i) { 461b1a0a8a3SJed Brown ierr = KSPSolve(osm->ksp[i],osm->x[i],osm->y[i]);CHKERRQ(ierr); 462b1a0a8a3SJed Brown } 463f746d493SDmitry Karpeev ierr = VecScatterBegin(osm->gprolongation,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); 464f746d493SDmitry Karpeev ierr = VecScatterEnd(osm->gprolongation,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); 465b1a0a8a3SJed Brown PetscFunctionReturn(0); 466b1a0a8a3SJed Brown } 467b1a0a8a3SJed Brown 468b1a0a8a3SJed Brown #undef __FUNCT__ 469f746d493SDmitry Karpeev #define __FUNCT__ "PCApplyTranspose_GASM" 470f746d493SDmitry Karpeev static PetscErrorCode PCApplyTranspose_GASM(PC pc,Vec x,Vec y) 471b1a0a8a3SJed Brown { 472f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 473b1a0a8a3SJed Brown PetscErrorCode ierr; 474f746d493SDmitry Karpeev PetscInt i; 475b1a0a8a3SJed Brown ScatterMode forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE; 476b1a0a8a3SJed Brown 477b1a0a8a3SJed Brown PetscFunctionBegin; 478b1a0a8a3SJed Brown /* 479b1a0a8a3SJed Brown Support for limiting the restriction or interpolation to only local 480b1a0a8a3SJed Brown subdomain values (leaving the other values 0). 481b1a0a8a3SJed Brown 482f746d493SDmitry Karpeev Note: these are reversed from the PCApply_GASM() because we are applying the 483b1a0a8a3SJed Brown transpose of the three terms 484b1a0a8a3SJed Brown */ 485f746d493SDmitry Karpeev if (!(osm->type & PC_GASM_INTERPOLATE)) { 486b1a0a8a3SJed Brown forward = SCATTER_FORWARD_LOCAL; 487b1a0a8a3SJed Brown /* have to zero the work RHS since scatter may leave some slots empty */ 488f746d493SDmitry Karpeev ierr = VecZeroEntries(osm->gx);CHKERRQ(ierr); 489b1a0a8a3SJed Brown } 490f746d493SDmitry Karpeev if (!(osm->type & PC_GASM_RESTRICT)) { 491b1a0a8a3SJed Brown reverse = SCATTER_REVERSE_LOCAL; 492b1a0a8a3SJed Brown } 493b1a0a8a3SJed Brown 494ab3e923fSDmitry Karpeev ierr = VecScatterBegin(osm->grestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 495b1a0a8a3SJed Brown ierr = VecZeroEntries(y);CHKERRQ(ierr); 496ab3e923fSDmitry Karpeev ierr = VecScatterEnd(osm->grestriction,x,osm->gx,INSERT_VALUES,forward);CHKERRQ(ierr); 497b1a0a8a3SJed Brown /* do the local solves */ 498ab3e923fSDmitry 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. */ 499b1a0a8a3SJed Brown ierr = KSPSolveTranspose(osm->ksp[i],osm->x[i],osm->y[i]);CHKERRQ(ierr); 500b1a0a8a3SJed Brown } 501ab3e923fSDmitry Karpeev ierr = VecScatterBegin(osm->gprolongation,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); 502ab3e923fSDmitry Karpeev ierr = VecScatterEnd(osm->gprolongation,osm->gy,y,ADD_VALUES,reverse);CHKERRQ(ierr); 503b1a0a8a3SJed Brown PetscFunctionReturn(0); 504b1a0a8a3SJed Brown } 505b1a0a8a3SJed Brown 506b1a0a8a3SJed Brown #undef __FUNCT__ 507a06653b4SBarry Smith #define __FUNCT__ "PCReset_GASM" 508a06653b4SBarry Smith static PetscErrorCode PCReset_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; 515b1a0a8a3SJed Brown if (osm->ksp) { 516f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 517a06653b4SBarry Smith ierr = KSPReset(osm->ksp[i]);CHKERRQ(ierr); 518b1a0a8a3SJed Brown } 519b1a0a8a3SJed Brown } 520b1a0a8a3SJed Brown if (osm->pmat) { 521f746d493SDmitry Karpeev if (osm->n > 0) { 522ab3e923fSDmitry Karpeev ierr = MatDestroyMatrices(osm->n,&osm->pmat);CHKERRQ(ierr); 523b1a0a8a3SJed Brown } 524b1a0a8a3SJed Brown } 525d34fcf5fSBarry Smith if (osm->x) { 526f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 527d34fcf5fSBarry Smith if (osm->x[i]) {ierr = VecDestroy(osm->x[i]);CHKERRQ(ierr);} 528d34fcf5fSBarry Smith if (osm->y[i]) {ierr = VecDestroy(osm->y[i]);CHKERRQ(ierr);} 529b1a0a8a3SJed Brown } 530d34fcf5fSBarry Smith } 531d34fcf5fSBarry Smith if (osm->gx) {ierr = VecDestroy(osm->gx);CHKERRQ(ierr);} 532d34fcf5fSBarry Smith if (osm->gy) {ierr = VecDestroy(osm->gy);CHKERRQ(ierr);} 533ab3e923fSDmitry Karpeev 534d34fcf5fSBarry Smith if (osm->grestriction) {ierr = VecScatterDestroy(osm->grestriction);CHKERRQ(ierr);} 535d34fcf5fSBarry Smith if (osm->gprolongation) {ierr = VecScatterDestroy(osm->gprolongation);CHKERRQ(ierr);} 536d34fcf5fSBarry Smith if (osm->is) {ierr = PCGASMDestroySubdomains(osm->n,osm->is,osm->is_local);CHKERRQ(ierr); osm->is = 0;} 537d34fcf5fSBarry Smith if (osm->gis) {ierr = ISDestroy(osm->gis);CHKERRQ(ierr);} 538d34fcf5fSBarry Smith if (osm->gis_local) {ierr = ISDestroy(osm->gis_local);CHKERRQ(ierr);} 539a06653b4SBarry Smith PetscFunctionReturn(0); 540a06653b4SBarry Smith } 541a06653b4SBarry Smith 542a06653b4SBarry Smith #undef __FUNCT__ 543a06653b4SBarry Smith #define __FUNCT__ "PCDestroy_GASM" 544a06653b4SBarry Smith static PetscErrorCode PCDestroy_GASM(PC pc) 545a06653b4SBarry Smith { 546a06653b4SBarry Smith PC_GASM *osm = (PC_GASM*)pc->data; 547a06653b4SBarry Smith PetscErrorCode ierr; 548a06653b4SBarry Smith PetscInt i; 549a06653b4SBarry Smith 550a06653b4SBarry Smith PetscFunctionBegin; 551a06653b4SBarry Smith ierr = PCReset_GASM(pc);CHKERRQ(ierr); 552a06653b4SBarry Smith if (osm->ksp) { 553a06653b4SBarry Smith for (i=0; i<osm->n; i++) { 554a06653b4SBarry Smith ierr = KSPDestroy(osm->ksp[i]);CHKERRQ(ierr); 555a06653b4SBarry Smith } 556a06653b4SBarry Smith ierr = PetscFree(osm->ksp);CHKERRQ(ierr); 557a06653b4SBarry Smith } 558a06653b4SBarry Smith ierr = PetscFree(osm->x);CHKERRQ(ierr); 559a06653b4SBarry Smith ierr = PetscFree(osm->y);CHKERRQ(ierr); 560c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 561b1a0a8a3SJed Brown PetscFunctionReturn(0); 562b1a0a8a3SJed Brown } 563b1a0a8a3SJed Brown 564b1a0a8a3SJed Brown #undef __FUNCT__ 565f746d493SDmitry Karpeev #define __FUNCT__ "PCSetFromOptions_GASM" 566ab3e923fSDmitry Karpeev static PetscErrorCode PCSetFromOptions_GASM(PC pc) { 567f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 568b1a0a8a3SJed Brown PetscErrorCode ierr; 569b1a0a8a3SJed Brown PetscInt blocks,ovl; 570b1a0a8a3SJed Brown PetscBool symset,flg; 571f746d493SDmitry Karpeev PCGASMType gasmtype; 572b1a0a8a3SJed Brown 573b1a0a8a3SJed Brown PetscFunctionBegin; 574b1a0a8a3SJed Brown /* set the type to symmetric if matrix is symmetric */ 575b1a0a8a3SJed Brown if (!osm->type_set && pc->pmat) { 576b1a0a8a3SJed Brown ierr = MatIsSymmetricKnown(pc->pmat,&symset,&flg);CHKERRQ(ierr); 577f746d493SDmitry Karpeev if (symset && flg) { osm->type = PC_GASM_BASIC; } 578b1a0a8a3SJed Brown } 57944fe18e5SDmitry Karpeev ierr = PetscOptionsHead("Generalized additive Schwarz options");CHKERRQ(ierr); 580f746d493SDmitry Karpeev ierr = PetscOptionsInt("-pc_gasm_blocks","Number of subdomains","PCGASMSetTotalSubdomains",osm->n,&blocks,&flg);CHKERRQ(ierr); 581ab3e923fSDmitry Karpeev if (flg) {ierr = PCGASMSetTotalSubdomains(pc,blocks);CHKERRQ(ierr); } 582f746d493SDmitry Karpeev ierr = PetscOptionsInt("-pc_gasm_overlap","Number of grid points overlap","PCGASMSetOverlap",osm->overlap,&ovl,&flg);CHKERRQ(ierr); 583f746d493SDmitry Karpeev if (flg) {ierr = PCGASMSetOverlap(pc,ovl);CHKERRQ(ierr); } 584b1a0a8a3SJed Brown flg = PETSC_FALSE; 585f746d493SDmitry Karpeev ierr = PetscOptionsEnum("-pc_gasm_type","Type of restriction/extension","PCGASMSetType",PCGASMTypes,(PetscEnum)osm->type,(PetscEnum*)&gasmtype,&flg);CHKERRQ(ierr); 586f746d493SDmitry Karpeev if (flg) {ierr = PCGASMSetType(pc,gasmtype);CHKERRQ(ierr); } 587b1a0a8a3SJed Brown ierr = PetscOptionsTail();CHKERRQ(ierr); 588b1a0a8a3SJed Brown PetscFunctionReturn(0); 589b1a0a8a3SJed Brown } 590b1a0a8a3SJed Brown 591b1a0a8a3SJed Brown /*------------------------------------------------------------------------------------*/ 592b1a0a8a3SJed Brown 593b1a0a8a3SJed Brown EXTERN_C_BEGIN 594b1a0a8a3SJed Brown #undef __FUNCT__ 595f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetLocalSubdomains_GASM" 5967087cfbeSBarry Smith PetscErrorCode PCGASMSetLocalSubdomains_GASM(PC pc,PetscInt n,IS is[],IS is_local[]) 597b1a0a8a3SJed Brown { 598f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 599b1a0a8a3SJed Brown PetscErrorCode ierr; 600b1a0a8a3SJed Brown PetscInt i; 601b1a0a8a3SJed Brown 602b1a0a8a3SJed Brown PetscFunctionBegin; 603b1a0a8a3SJed Brown if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Each process must have 1 or more blocks, n = %D",n); 604f746d493SDmitry Karpeev if (pc->setupcalled && (n != osm->n || is)) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetLocalSubdomains() should be called before calling PCSetUp()."); 605b1a0a8a3SJed Brown 606b1a0a8a3SJed Brown if (!pc->setupcalled) { 60793c1ec32SDmitry Karpeev osm->n = n; 60893c1ec32SDmitry Karpeev osm->is = 0; 60993c1ec32SDmitry Karpeev osm->is_local = 0; 610b1a0a8a3SJed Brown if (is) { 611b1a0a8a3SJed Brown for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)is[i]);CHKERRQ(ierr);} 612b1a0a8a3SJed Brown } 613b1a0a8a3SJed Brown if (is_local) { 614b1a0a8a3SJed Brown for (i=0; i<n; i++) {ierr = PetscObjectReference((PetscObject)is_local[i]);CHKERRQ(ierr);} 615b1a0a8a3SJed Brown } 616b1a0a8a3SJed Brown if (osm->is) { 617ab3e923fSDmitry Karpeev ierr = PCGASMDestroySubdomains(osm->n,osm->is,osm->is_local);CHKERRQ(ierr); 618b1a0a8a3SJed Brown } 619b1a0a8a3SJed Brown if (is) { 620b1a0a8a3SJed Brown ierr = PetscMalloc(n*sizeof(IS),&osm->is);CHKERRQ(ierr); 621b1a0a8a3SJed Brown for (i=0; i<n; i++) { osm->is[i] = is[i]; } 622f746d493SDmitry Karpeev /* Flag indicating that the user has set overlapping subdomains so PCGASM should not increase their size. */ 623b1a0a8a3SJed Brown osm->overlap = -1; 624b1a0a8a3SJed Brown } 625b1a0a8a3SJed Brown if (is_local) { 626b1a0a8a3SJed Brown ierr = PetscMalloc(n*sizeof(IS),&osm->is_local);CHKERRQ(ierr); 627b1a0a8a3SJed Brown for (i=0; i<n; i++) { osm->is_local[i] = is_local[i]; } 628b1a0a8a3SJed Brown } 629b1a0a8a3SJed Brown } 630b1a0a8a3SJed Brown PetscFunctionReturn(0); 631b1a0a8a3SJed Brown } 632b1a0a8a3SJed Brown EXTERN_C_END 633b1a0a8a3SJed Brown 634b1a0a8a3SJed Brown EXTERN_C_BEGIN 635b1a0a8a3SJed Brown #undef __FUNCT__ 636f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetTotalSubdomains_GASM" 6377087cfbeSBarry Smith PetscErrorCode PCGASMSetTotalSubdomains_GASM(PC pc,PetscInt N) { 638f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 639b1a0a8a3SJed Brown PetscErrorCode ierr; 640b1a0a8a3SJed Brown PetscMPIInt rank,size; 641b1a0a8a3SJed Brown PetscInt n; 642b1a0a8a3SJed Brown 643b1a0a8a3SJed Brown PetscFunctionBegin; 644b1a0a8a3SJed Brown if (N < 1) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Number of total blocks must be > 0, N = %D",N); 645b1a0a8a3SJed Brown 646b1a0a8a3SJed Brown /* 647b1a0a8a3SJed Brown Split the subdomains equally among all processors 648b1a0a8a3SJed Brown */ 649b1a0a8a3SJed Brown ierr = MPI_Comm_rank(((PetscObject)pc)->comm,&rank);CHKERRQ(ierr); 650b1a0a8a3SJed Brown ierr = MPI_Comm_size(((PetscObject)pc)->comm,&size);CHKERRQ(ierr); 651b1a0a8a3SJed Brown n = N/size + ((N % size) > rank); 652b1a0a8a3SJed Brown if (!n) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Process %d must have at least one block: total processors %d total blocks %D",(int)rank,(int)size,N); 653f746d493SDmitry Karpeev if (pc->setupcalled && n != osm->n) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetTotalSubdomains() should be called before PCSetUp()."); 654b1a0a8a3SJed Brown if (!pc->setupcalled) { 655b1a0a8a3SJed Brown if (osm->is) { 656ab3e923fSDmitry Karpeev ierr = PCGASMDestroySubdomains(osm->n,osm->is,osm->is_local);CHKERRQ(ierr); 657b1a0a8a3SJed Brown } 65893c1ec32SDmitry Karpeev osm->N = N; 659f746d493SDmitry Karpeev osm->n = n; 660b1a0a8a3SJed Brown osm->is = 0; 661b1a0a8a3SJed Brown osm->is_local = 0; 662b1a0a8a3SJed Brown } 663b1a0a8a3SJed Brown PetscFunctionReturn(0); 664ab3e923fSDmitry Karpeev }/* PCGASMSetTotalSubdomains_GASM() */ 665b1a0a8a3SJed Brown EXTERN_C_END 666b1a0a8a3SJed Brown 667b1a0a8a3SJed Brown EXTERN_C_BEGIN 668b1a0a8a3SJed Brown #undef __FUNCT__ 669f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetOverlap_GASM" 6707087cfbeSBarry Smith PetscErrorCode PCGASMSetOverlap_GASM(PC pc,PetscInt ovl) 671b1a0a8a3SJed Brown { 672f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 673b1a0a8a3SJed Brown 674b1a0a8a3SJed Brown PetscFunctionBegin; 675b1a0a8a3SJed Brown if (ovl < 0) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap value requested"); 676f746d493SDmitry Karpeev if (pc->setupcalled && ovl != osm->overlap) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetOverlap() should be called before PCSetUp()."); 677b1a0a8a3SJed Brown if (!pc->setupcalled) { 678b1a0a8a3SJed Brown osm->overlap = ovl; 679b1a0a8a3SJed Brown } 680b1a0a8a3SJed Brown PetscFunctionReturn(0); 681b1a0a8a3SJed Brown } 682b1a0a8a3SJed Brown EXTERN_C_END 683b1a0a8a3SJed Brown 684b1a0a8a3SJed Brown EXTERN_C_BEGIN 685b1a0a8a3SJed Brown #undef __FUNCT__ 686f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetType_GASM" 6877087cfbeSBarry Smith PetscErrorCode PCGASMSetType_GASM(PC pc,PCGASMType type) 688b1a0a8a3SJed Brown { 689f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 690b1a0a8a3SJed Brown 691b1a0a8a3SJed Brown PetscFunctionBegin; 692b1a0a8a3SJed Brown osm->type = type; 693b1a0a8a3SJed Brown osm->type_set = PETSC_TRUE; 694b1a0a8a3SJed Brown PetscFunctionReturn(0); 695b1a0a8a3SJed Brown } 696b1a0a8a3SJed Brown EXTERN_C_END 697b1a0a8a3SJed Brown 698b1a0a8a3SJed Brown EXTERN_C_BEGIN 699b1a0a8a3SJed Brown #undef __FUNCT__ 700f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetSortIndices_GASM" 7017087cfbeSBarry Smith PetscErrorCode PCGASMSetSortIndices_GASM(PC pc,PetscBool doSort) 702b1a0a8a3SJed Brown { 703f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 704b1a0a8a3SJed Brown 705b1a0a8a3SJed Brown PetscFunctionBegin; 706b1a0a8a3SJed Brown osm->sort_indices = doSort; 707b1a0a8a3SJed Brown PetscFunctionReturn(0); 708b1a0a8a3SJed Brown } 709b1a0a8a3SJed Brown EXTERN_C_END 710b1a0a8a3SJed Brown 711b1a0a8a3SJed Brown EXTERN_C_BEGIN 712b1a0a8a3SJed Brown #undef __FUNCT__ 713f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMGetSubKSP_GASM" 71444fe18e5SDmitry Karpeev /* 71544fe18e5SDmitry Karpeev FIX: This routine might need to be modified once multiple ranks per subdomain are allowed. 71644fe18e5SDmitry Karpeev In particular, it would upset the global subdomain number calculation. 71744fe18e5SDmitry Karpeev */ 7187087cfbeSBarry Smith PetscErrorCode PCGASMGetSubKSP_GASM(PC pc,PetscInt *n,PetscInt *first,KSP **ksp) 719b1a0a8a3SJed Brown { 720f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 721b1a0a8a3SJed Brown PetscErrorCode ierr; 722b1a0a8a3SJed Brown 723b1a0a8a3SJed Brown PetscFunctionBegin; 724ab3e923fSDmitry 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"); 725b1a0a8a3SJed Brown 726ab3e923fSDmitry Karpeev if (n) { 727ab3e923fSDmitry Karpeev *n = osm->n; 728b1a0a8a3SJed Brown } 729ab3e923fSDmitry Karpeev if (first) { 730ab3e923fSDmitry Karpeev ierr = MPI_Scan(&osm->n,first,1,MPIU_INT,MPI_SUM,((PetscObject)pc)->comm);CHKERRQ(ierr); 731ab3e923fSDmitry Karpeev *first -= osm->n; 732b1a0a8a3SJed Brown } 733b1a0a8a3SJed Brown if (ksp) { 734b1a0a8a3SJed Brown /* Assume that local solves are now different; not necessarily 735f746d493SDmitry Karpeev true though! This flag is used only for PCView_GASM() */ 736b1a0a8a3SJed Brown *ksp = osm->ksp; 737b1a0a8a3SJed Brown osm->same_local_solves = PETSC_FALSE; 738b1a0a8a3SJed Brown } 739b1a0a8a3SJed Brown PetscFunctionReturn(0); 740ab3e923fSDmitry Karpeev }/* PCGASMGetSubKSP_GASM() */ 741b1a0a8a3SJed Brown EXTERN_C_END 742b1a0a8a3SJed Brown 743b1a0a8a3SJed Brown 744b1a0a8a3SJed Brown #undef __FUNCT__ 745f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetLocalSubdomains" 746b1a0a8a3SJed Brown /*@C 747f746d493SDmitry Karpeev PCGASMSetLocalSubdomains - Sets the local subdomains (for this processor 748b1a0a8a3SJed Brown only) for the additive Schwarz preconditioner. 749b1a0a8a3SJed Brown 750b1a0a8a3SJed Brown Collective on PC 751b1a0a8a3SJed Brown 752b1a0a8a3SJed Brown Input Parameters: 753b1a0a8a3SJed Brown + pc - the preconditioner context 754b1a0a8a3SJed Brown . n - the number of subdomains for this processor (default value = 1) 755b1a0a8a3SJed Brown . is - the index set that defines the subdomains for this processor 756b1a0a8a3SJed Brown (or PETSC_NULL for PETSc to determine subdomains) 757b1a0a8a3SJed Brown - is_local - the index sets that define the local part of the subdomains for this processor 758b1a0a8a3SJed Brown (or PETSC_NULL to use the default of 1 subdomain per process) 759b1a0a8a3SJed Brown 760b1a0a8a3SJed Brown Notes: 761b1a0a8a3SJed Brown The IS numbering is in the parallel, global numbering of the vector. 762b1a0a8a3SJed Brown 763f746d493SDmitry Karpeev By default the GASM preconditioner uses 1 block per processor. 764b1a0a8a3SJed Brown 765f746d493SDmitry Karpeev Use PCGASMSetTotalSubdomains() to set the subdomains for all processors. 766b1a0a8a3SJed Brown 767b1a0a8a3SJed Brown Level: advanced 768b1a0a8a3SJed Brown 769f746d493SDmitry Karpeev .keywords: PC, GASM, set, local, subdomains, additive Schwarz 770b1a0a8a3SJed Brown 771f746d493SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetOverlap(), PCGASMGetSubKSP(), 772f746d493SDmitry Karpeev PCGASMCreateSubdomains2D(), PCGASMGetLocalSubdomains() 773b1a0a8a3SJed Brown @*/ 7747087cfbeSBarry Smith PetscErrorCode PCGASMSetLocalSubdomains(PC pc,PetscInt n,IS is[],IS is_local[]) 775b1a0a8a3SJed Brown { 776b1a0a8a3SJed Brown PetscErrorCode ierr; 777b1a0a8a3SJed Brown 778b1a0a8a3SJed Brown PetscFunctionBegin; 779b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 780f746d493SDmitry Karpeev ierr = PetscTryMethod(pc,"PCGASMSetLocalSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,n,is,is_local));CHKERRQ(ierr); 781b1a0a8a3SJed Brown PetscFunctionReturn(0); 782b1a0a8a3SJed Brown } 783b1a0a8a3SJed Brown 784b1a0a8a3SJed Brown #undef __FUNCT__ 785f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetTotalSubdomains" 786b1a0a8a3SJed Brown /*@C 787f746d493SDmitry Karpeev PCGASMSetTotalSubdomains - Sets the subdomains for all processor for the 788b1a0a8a3SJed Brown additive Schwarz preconditioner. Either all or no processors in the 789b1a0a8a3SJed Brown PC communicator must call this routine, with the same index sets. 790b1a0a8a3SJed Brown 791b1a0a8a3SJed Brown Collective on PC 792b1a0a8a3SJed Brown 793b1a0a8a3SJed Brown Input Parameters: 794b1a0a8a3SJed Brown + pc - the preconditioner context 795b1a0a8a3SJed Brown . n - the number of subdomains for all processors 796b1a0a8a3SJed Brown . is - the index sets that define the subdomains for all processor 797b1a0a8a3SJed Brown (or PETSC_NULL for PETSc to determine subdomains) 798b1a0a8a3SJed Brown - is_local - the index sets that define the local part of the subdomains for this processor 799b1a0a8a3SJed Brown (or PETSC_NULL to use the default of 1 subdomain per process) 800b1a0a8a3SJed Brown 801b1a0a8a3SJed Brown Options Database Key: 802b1a0a8a3SJed Brown To set the total number of subdomain blocks rather than specify the 803b1a0a8a3SJed Brown index sets, use the option 804f746d493SDmitry Karpeev . -pc_gasm_blocks <blks> - Sets total blocks 805b1a0a8a3SJed Brown 806b1a0a8a3SJed Brown Notes: 807b1a0a8a3SJed Brown Currently you cannot use this to set the actual subdomains with the argument is. 808b1a0a8a3SJed Brown 809f746d493SDmitry Karpeev By default the GASM preconditioner uses 1 block per processor. 810b1a0a8a3SJed Brown 811b1a0a8a3SJed Brown These index sets cannot be destroyed until after completion of the 812f746d493SDmitry Karpeev linear solves for which the GASM preconditioner is being used. 813b1a0a8a3SJed Brown 814f746d493SDmitry Karpeev Use PCGASMSetLocalSubdomains() to set local subdomains. 815b1a0a8a3SJed Brown 816b1a0a8a3SJed Brown Level: advanced 817b1a0a8a3SJed Brown 818f746d493SDmitry Karpeev .keywords: PC, GASM, set, total, global, subdomains, additive Schwarz 819b1a0a8a3SJed Brown 820f746d493SDmitry Karpeev .seealso: PCGASMSetLocalSubdomains(), PCGASMSetOverlap(), PCGASMGetSubKSP(), 821f746d493SDmitry Karpeev PCGASMCreateSubdomains2D() 822b1a0a8a3SJed Brown @*/ 8237087cfbeSBarry Smith PetscErrorCode PCGASMSetTotalSubdomains(PC pc,PetscInt N) 824b1a0a8a3SJed Brown { 825b1a0a8a3SJed Brown PetscErrorCode ierr; 826b1a0a8a3SJed Brown 827b1a0a8a3SJed Brown PetscFunctionBegin; 828b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 829f746d493SDmitry Karpeev ierr = PetscTryMethod(pc,"PCGASMSetTotalSubdomains_C",(PC,PetscInt),(pc,N));CHKERRQ(ierr); 830b1a0a8a3SJed Brown PetscFunctionReturn(0); 831b1a0a8a3SJed Brown } 832b1a0a8a3SJed Brown 833b1a0a8a3SJed Brown #undef __FUNCT__ 834f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetOverlap" 835b1a0a8a3SJed Brown /*@ 836f746d493SDmitry Karpeev PCGASMSetOverlap - Sets the overlap between a pair of subdomains for the 837b1a0a8a3SJed Brown additive Schwarz preconditioner. Either all or no processors in the 838b1a0a8a3SJed Brown PC communicator must call this routine. 839b1a0a8a3SJed Brown 840b1a0a8a3SJed Brown Logically Collective on PC 841b1a0a8a3SJed Brown 842b1a0a8a3SJed Brown Input Parameters: 843b1a0a8a3SJed Brown + pc - the preconditioner context 844b1a0a8a3SJed Brown - ovl - the amount of overlap between subdomains (ovl >= 0, default value = 1) 845b1a0a8a3SJed Brown 846b1a0a8a3SJed Brown Options Database Key: 847f746d493SDmitry Karpeev . -pc_gasm_overlap <ovl> - Sets overlap 848b1a0a8a3SJed Brown 849b1a0a8a3SJed Brown Notes: 850f746d493SDmitry Karpeev By default the GASM preconditioner uses 1 block per processor. To use 851f746d493SDmitry Karpeev multiple blocks per perocessor, see PCGASMSetTotalSubdomains() and 852f746d493SDmitry Karpeev PCGASMSetLocalSubdomains() (and the option -pc_gasm_blocks <blks>). 853b1a0a8a3SJed Brown 854b1a0a8a3SJed Brown The overlap defaults to 1, so if one desires that no additional 855b1a0a8a3SJed Brown overlap be computed beyond what may have been set with a call to 856f746d493SDmitry Karpeev PCGASMSetTotalSubdomains() or PCGASMSetLocalSubdomains(), then ovl 857b1a0a8a3SJed Brown must be set to be 0. In particular, if one does not explicitly set 858b1a0a8a3SJed Brown the subdomains an application code, then all overlap would be computed 859f746d493SDmitry Karpeev internally by PETSc, and using an overlap of 0 would result in an GASM 860b1a0a8a3SJed Brown variant that is equivalent to the block Jacobi preconditioner. 861b1a0a8a3SJed Brown 862b1a0a8a3SJed Brown Note that one can define initial index sets with any overlap via 863f746d493SDmitry Karpeev PCGASMSetTotalSubdomains() or PCGASMSetLocalSubdomains(); the routine 864f746d493SDmitry Karpeev PCGASMSetOverlap() merely allows PETSc to extend that overlap further 865b1a0a8a3SJed Brown if desired. 866b1a0a8a3SJed Brown 867b1a0a8a3SJed Brown Level: intermediate 868b1a0a8a3SJed Brown 869f746d493SDmitry Karpeev .keywords: PC, GASM, set, overlap 870b1a0a8a3SJed Brown 871f746d493SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetLocalSubdomains(), PCGASMGetSubKSP(), 872f746d493SDmitry Karpeev PCGASMCreateSubdomains2D(), PCGASMGetLocalSubdomains() 873b1a0a8a3SJed Brown @*/ 8747087cfbeSBarry Smith PetscErrorCode PCGASMSetOverlap(PC pc,PetscInt ovl) 875b1a0a8a3SJed Brown { 876b1a0a8a3SJed Brown PetscErrorCode ierr; 877b1a0a8a3SJed Brown 878b1a0a8a3SJed Brown PetscFunctionBegin; 879b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 880b1a0a8a3SJed Brown PetscValidLogicalCollectiveInt(pc,ovl,2); 881f746d493SDmitry Karpeev ierr = PetscTryMethod(pc,"PCGASMSetOverlap_C",(PC,PetscInt),(pc,ovl));CHKERRQ(ierr); 882b1a0a8a3SJed Brown PetscFunctionReturn(0); 883b1a0a8a3SJed Brown } 884b1a0a8a3SJed Brown 885b1a0a8a3SJed Brown #undef __FUNCT__ 886f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetType" 887b1a0a8a3SJed Brown /*@ 888f746d493SDmitry Karpeev PCGASMSetType - Sets the type of restriction and interpolation used 889b1a0a8a3SJed Brown for local problems in the additive Schwarz method. 890b1a0a8a3SJed Brown 891b1a0a8a3SJed Brown Logically Collective on PC 892b1a0a8a3SJed Brown 893b1a0a8a3SJed Brown Input Parameters: 894b1a0a8a3SJed Brown + pc - the preconditioner context 895f746d493SDmitry Karpeev - type - variant of GASM, one of 896b1a0a8a3SJed Brown .vb 897f746d493SDmitry Karpeev PC_GASM_BASIC - full interpolation and restriction 898f746d493SDmitry Karpeev PC_GASM_RESTRICT - full restriction, local processor interpolation 899f746d493SDmitry Karpeev PC_GASM_INTERPOLATE - full interpolation, local processor restriction 900f746d493SDmitry Karpeev PC_GASM_NONE - local processor restriction and interpolation 901b1a0a8a3SJed Brown .ve 902b1a0a8a3SJed Brown 903b1a0a8a3SJed Brown Options Database Key: 904f746d493SDmitry Karpeev . -pc_gasm_type [basic,restrict,interpolate,none] - Sets GASM type 905b1a0a8a3SJed Brown 906b1a0a8a3SJed Brown Level: intermediate 907b1a0a8a3SJed Brown 908f746d493SDmitry Karpeev .keywords: PC, GASM, set, type 909b1a0a8a3SJed Brown 910f746d493SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetTotalSubdomains(), PCGASMGetSubKSP(), 911f746d493SDmitry Karpeev PCGASMCreateSubdomains2D() 912b1a0a8a3SJed Brown @*/ 9137087cfbeSBarry Smith PetscErrorCode PCGASMSetType(PC pc,PCGASMType type) 914b1a0a8a3SJed Brown { 915b1a0a8a3SJed Brown PetscErrorCode ierr; 916b1a0a8a3SJed Brown 917b1a0a8a3SJed Brown PetscFunctionBegin; 918b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 919b1a0a8a3SJed Brown PetscValidLogicalCollectiveEnum(pc,type,2); 920f746d493SDmitry Karpeev ierr = PetscTryMethod(pc,"PCGASMSetType_C",(PC,PCGASMType),(pc,type));CHKERRQ(ierr); 921b1a0a8a3SJed Brown PetscFunctionReturn(0); 922b1a0a8a3SJed Brown } 923b1a0a8a3SJed Brown 924b1a0a8a3SJed Brown #undef __FUNCT__ 925f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMSetSortIndices" 926b1a0a8a3SJed Brown /*@ 927f746d493SDmitry Karpeev PCGASMSetSortIndices - Determines whether subdomain indices are sorted. 928b1a0a8a3SJed Brown 929b1a0a8a3SJed Brown Logically Collective on PC 930b1a0a8a3SJed Brown 931b1a0a8a3SJed Brown Input Parameters: 932b1a0a8a3SJed Brown + pc - the preconditioner context 933b1a0a8a3SJed Brown - doSort - sort the subdomain indices 934b1a0a8a3SJed Brown 935b1a0a8a3SJed Brown Level: intermediate 936b1a0a8a3SJed Brown 937f746d493SDmitry Karpeev .keywords: PC, GASM, set, type 938b1a0a8a3SJed Brown 939f746d493SDmitry Karpeev .seealso: PCGASMSetLocalSubdomains(), PCGASMSetTotalSubdomains(), PCGASMGetSubKSP(), 940f746d493SDmitry Karpeev PCGASMCreateSubdomains2D() 941b1a0a8a3SJed Brown @*/ 9427087cfbeSBarry Smith PetscErrorCode PCGASMSetSortIndices(PC pc,PetscBool doSort) 943b1a0a8a3SJed Brown { 944b1a0a8a3SJed Brown PetscErrorCode ierr; 945b1a0a8a3SJed Brown 946b1a0a8a3SJed Brown PetscFunctionBegin; 947b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 948b1a0a8a3SJed Brown PetscValidLogicalCollectiveBool(pc,doSort,2); 949f746d493SDmitry Karpeev ierr = PetscTryMethod(pc,"PCGASMSetSortIndices_C",(PC,PetscBool),(pc,doSort));CHKERRQ(ierr); 950b1a0a8a3SJed Brown PetscFunctionReturn(0); 951b1a0a8a3SJed Brown } 952b1a0a8a3SJed Brown 953b1a0a8a3SJed Brown #undef __FUNCT__ 954f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMGetSubKSP" 955b1a0a8a3SJed Brown /*@C 956f746d493SDmitry Karpeev PCGASMGetSubKSP - Gets the local KSP contexts for all blocks on 957b1a0a8a3SJed Brown this processor. 958b1a0a8a3SJed Brown 959b1a0a8a3SJed Brown Collective on PC iff first_local is requested 960b1a0a8a3SJed Brown 961b1a0a8a3SJed Brown Input Parameter: 962b1a0a8a3SJed Brown . pc - the preconditioner context 963b1a0a8a3SJed Brown 964b1a0a8a3SJed Brown Output Parameters: 965b1a0a8a3SJed Brown + n_local - the number of blocks on this processor or PETSC_NULL 966b1a0a8a3SJed Brown . first_local - the global number of the first block on this processor or PETSC_NULL, 967b1a0a8a3SJed Brown all processors must request or all must pass PETSC_NULL 968b1a0a8a3SJed Brown - ksp - the array of KSP contexts 969b1a0a8a3SJed Brown 970b1a0a8a3SJed Brown Note: 971f746d493SDmitry Karpeev After PCGASMGetSubKSP() the array of KSPes is not to be freed 972b1a0a8a3SJed Brown 973b1a0a8a3SJed Brown Currently for some matrix implementations only 1 block per processor 974b1a0a8a3SJed Brown is supported. 975b1a0a8a3SJed Brown 976f746d493SDmitry Karpeev You must call KSPSetUp() before calling PCGASMGetSubKSP(). 977b1a0a8a3SJed Brown 978b1a0a8a3SJed Brown Level: advanced 979b1a0a8a3SJed Brown 980f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, get, sub, KSP, context 981b1a0a8a3SJed Brown 982f746d493SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetTotalSubdomains(), PCGASMSetOverlap(), 983f746d493SDmitry Karpeev PCGASMCreateSubdomains2D(), 984b1a0a8a3SJed Brown @*/ 9857087cfbeSBarry Smith PetscErrorCode PCGASMGetSubKSP(PC pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp[]) 986b1a0a8a3SJed Brown { 987b1a0a8a3SJed Brown PetscErrorCode ierr; 988b1a0a8a3SJed Brown 989b1a0a8a3SJed Brown PetscFunctionBegin; 990b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 991f746d493SDmitry Karpeev ierr = PetscUseMethod(pc,"PCGASMGetSubKSP_C",(PC,PetscInt*,PetscInt*,KSP **),(pc,n_local,first_local,ksp));CHKERRQ(ierr); 992b1a0a8a3SJed Brown PetscFunctionReturn(0); 993b1a0a8a3SJed Brown } 994b1a0a8a3SJed Brown 995b1a0a8a3SJed Brown /* -------------------------------------------------------------------------------------*/ 996b1a0a8a3SJed Brown /*MC 997f746d493SDmitry Karpeev PCGASM - Use the (restricted) additive Schwarz method, each block is (approximately) solved with 998b1a0a8a3SJed Brown its own KSP object. 999b1a0a8a3SJed Brown 1000b1a0a8a3SJed Brown Options Database Keys: 1001f746d493SDmitry Karpeev + -pc_gasm_truelocal - Activates PCGGASMSetUseTrueLocal() 1002f746d493SDmitry Karpeev . -pc_gasm_blocks <blks> - Sets total blocks 1003f746d493SDmitry Karpeev . -pc_gasm_overlap <ovl> - Sets overlap 1004f746d493SDmitry Karpeev - -pc_gasm_type [basic,restrict,interpolate,none] - Sets GASM type 1005b1a0a8a3SJed Brown 1006b1a0a8a3SJed Brown IMPORTANT: If you run with, for example, 3 blocks on 1 processor or 3 blocks on 3 processors you 1007f746d493SDmitry Karpeev will get a different convergence rate due to the default option of -pc_gasm_type restrict. Use 1008f746d493SDmitry Karpeev -pc_gasm_type basic to use the standard GASM. 1009b1a0a8a3SJed Brown 1010b1a0a8a3SJed Brown Notes: Each processor can have one or more blocks, but a block cannot be shared by more 1011b1a0a8a3SJed Brown than one processor. Defaults to one block per processor. 1012b1a0a8a3SJed Brown 1013b1a0a8a3SJed Brown To set options on the solvers for each block append -sub_ to all the KSP, and PC 1014b1a0a8a3SJed Brown options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly 1015b1a0a8a3SJed Brown 1016f746d493SDmitry Karpeev To set the options on the solvers separate for each block call PCGASMGetSubKSP() 1017b1a0a8a3SJed Brown and set the options directly on the resulting KSP object (you can access its PC 1018b1a0a8a3SJed Brown with KSPGetPC()) 1019b1a0a8a3SJed Brown 1020b1a0a8a3SJed Brown 1021b1a0a8a3SJed Brown Level: beginner 1022b1a0a8a3SJed Brown 1023b1a0a8a3SJed Brown Concepts: additive Schwarz method 1024b1a0a8a3SJed Brown 1025b1a0a8a3SJed Brown References: 1026b1a0a8a3SJed Brown An additive variant of the Schwarz alternating method for the case of many subregions 1027b1a0a8a3SJed Brown M Dryja, OB Widlund - Courant Institute, New York University Technical report 1028b1a0a8a3SJed Brown 1029b1a0a8a3SJed Brown Domain Decompositions: Parallel Multilevel Methods for Elliptic Partial Differential Equations, 1030b1a0a8a3SJed Brown Barry Smith, Petter Bjorstad, and William Gropp, Cambridge University Press, ISBN 0-521-49589-X. 1031b1a0a8a3SJed Brown 1032b1a0a8a3SJed Brown .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 1033f746d493SDmitry Karpeev PCBJACOBI, PCGASMSetUseTrueLocal(), PCGASMGetSubKSP(), PCGASMSetLocalSubdomains(), 1034f746d493SDmitry Karpeev PCGASMSetTotalSubdomains(), PCSetModifySubmatrices(), PCGASMSetOverlap(), PCGASMSetType() 1035b1a0a8a3SJed Brown 1036b1a0a8a3SJed Brown M*/ 1037b1a0a8a3SJed Brown 1038b1a0a8a3SJed Brown EXTERN_C_BEGIN 1039b1a0a8a3SJed Brown #undef __FUNCT__ 1040f746d493SDmitry Karpeev #define __FUNCT__ "PCCreate_GASM" 10417087cfbeSBarry Smith PetscErrorCode PCCreate_GASM(PC pc) 1042b1a0a8a3SJed Brown { 1043b1a0a8a3SJed Brown PetscErrorCode ierr; 1044f746d493SDmitry Karpeev PC_GASM *osm; 1045b1a0a8a3SJed Brown 1046b1a0a8a3SJed Brown PetscFunctionBegin; 1047f746d493SDmitry Karpeev ierr = PetscNewLog(pc,PC_GASM,&osm);CHKERRQ(ierr); 1048ab3e923fSDmitry Karpeev osm->N = PETSC_DECIDE; 1049ab3e923fSDmitry Karpeev osm->n = 0; 1050ab3e923fSDmitry Karpeev osm->nmax = 0; 1051b1a0a8a3SJed Brown osm->overlap = 1; 1052b1a0a8a3SJed Brown osm->ksp = 0; 1053ab3e923fSDmitry Karpeev osm->grestriction = 0; 1054ab3e923fSDmitry Karpeev osm->gprolongation = 0; 1055ab3e923fSDmitry Karpeev osm->gx = 0; 1056ab3e923fSDmitry Karpeev osm->gy = 0; 1057b1a0a8a3SJed Brown osm->x = 0; 1058b1a0a8a3SJed Brown osm->y = 0; 1059b1a0a8a3SJed Brown osm->is = 0; 1060b1a0a8a3SJed Brown osm->is_local = 0; 1061b1a0a8a3SJed Brown osm->pmat = 0; 1062f746d493SDmitry Karpeev osm->type = PC_GASM_RESTRICT; 1063b1a0a8a3SJed Brown osm->same_local_solves = PETSC_TRUE; 1064b1a0a8a3SJed Brown osm->sort_indices = PETSC_TRUE; 1065b1a0a8a3SJed Brown 1066b1a0a8a3SJed Brown pc->data = (void*)osm; 1067f746d493SDmitry Karpeev pc->ops->apply = PCApply_GASM; 1068f746d493SDmitry Karpeev pc->ops->applytranspose = PCApplyTranspose_GASM; 1069f746d493SDmitry Karpeev pc->ops->setup = PCSetUp_GASM; 1070a06653b4SBarry Smith pc->ops->reset = PCReset_GASM; 1071f746d493SDmitry Karpeev pc->ops->destroy = PCDestroy_GASM; 1072f746d493SDmitry Karpeev pc->ops->setfromoptions = PCSetFromOptions_GASM; 1073f746d493SDmitry Karpeev pc->ops->setuponblocks = PCSetUpOnBlocks_GASM; 1074f746d493SDmitry Karpeev pc->ops->view = PCView_GASM; 1075b1a0a8a3SJed Brown pc->ops->applyrichardson = 0; 1076b1a0a8a3SJed Brown 1077f746d493SDmitry Karpeev ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMSetLocalSubdomains_C","PCGASMSetLocalSubdomains_GASM", 1078f746d493SDmitry Karpeev PCGASMSetLocalSubdomains_GASM);CHKERRQ(ierr); 1079f746d493SDmitry Karpeev ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMSetTotalSubdomains_C","PCGASMSetTotalSubdomains_GASM", 1080f746d493SDmitry Karpeev PCGASMSetTotalSubdomains_GASM);CHKERRQ(ierr); 1081f746d493SDmitry Karpeev ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMSetOverlap_C","PCGASMSetOverlap_GASM", 1082f746d493SDmitry Karpeev PCGASMSetOverlap_GASM);CHKERRQ(ierr); 1083f746d493SDmitry Karpeev ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMSetType_C","PCGASMSetType_GASM", 1084f746d493SDmitry Karpeev PCGASMSetType_GASM);CHKERRQ(ierr); 1085f746d493SDmitry Karpeev ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMSetSortIndices_C","PCGASMSetSortIndices_GASM", 1086f746d493SDmitry Karpeev PCGASMSetSortIndices_GASM);CHKERRQ(ierr); 1087f746d493SDmitry Karpeev ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCGASMGetSubKSP_C","PCGASMGetSubKSP_GASM", 1088f746d493SDmitry Karpeev PCGASMGetSubKSP_GASM);CHKERRQ(ierr); 1089b1a0a8a3SJed Brown PetscFunctionReturn(0); 1090b1a0a8a3SJed Brown } 1091b1a0a8a3SJed Brown EXTERN_C_END 1092b1a0a8a3SJed Brown 1093b1a0a8a3SJed Brown 1094b1a0a8a3SJed Brown #undef __FUNCT__ 1095f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMCreateSubdomains" 1096b1a0a8a3SJed Brown /*@C 1097f746d493SDmitry Karpeev PCGASMCreateSubdomains - Creates the index sets for the overlapping Schwarz 1098b1a0a8a3SJed Brown preconditioner for a any problem on a general grid. 1099b1a0a8a3SJed Brown 1100b1a0a8a3SJed Brown Collective 1101b1a0a8a3SJed Brown 1102b1a0a8a3SJed Brown Input Parameters: 1103b1a0a8a3SJed Brown + A - The global matrix operator 1104b1a0a8a3SJed Brown - n - the number of local blocks 1105b1a0a8a3SJed Brown 1106b1a0a8a3SJed Brown Output Parameters: 1107b1a0a8a3SJed Brown . outis - the array of index sets defining the subdomains 1108b1a0a8a3SJed Brown 1109b1a0a8a3SJed Brown Level: advanced 1110b1a0a8a3SJed Brown 1111f746d493SDmitry Karpeev Note: this generates nonoverlapping subdomains; the PCGASM will generate the overlap 1112f746d493SDmitry Karpeev from these if you use PCGASMSetLocalSubdomains() 1113b1a0a8a3SJed Brown 1114b1a0a8a3SJed Brown In the Fortran version you must provide the array outis[] already allocated of length n. 1115b1a0a8a3SJed Brown 1116f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, create, subdomains, unstructured grid 1117b1a0a8a3SJed Brown 1118f746d493SDmitry Karpeev .seealso: PCGASMSetLocalSubdomains(), PCGASMDestroySubdomains() 1119b1a0a8a3SJed Brown @*/ 11207087cfbeSBarry Smith PetscErrorCode PCGASMCreateSubdomains(Mat A, PetscInt n, IS* outis[]) 1121b1a0a8a3SJed Brown { 1122b1a0a8a3SJed Brown MatPartitioning mpart; 1123b1a0a8a3SJed Brown const char *prefix; 1124b1a0a8a3SJed Brown PetscErrorCode (*f)(Mat,PetscBool *,MatReuse,Mat*); 1125b1a0a8a3SJed Brown PetscMPIInt size; 1126b1a0a8a3SJed Brown PetscInt i,j,rstart,rend,bs; 1127b1a0a8a3SJed Brown PetscBool iscopy = PETSC_FALSE,isbaij = PETSC_FALSE,foundpart = PETSC_FALSE; 1128b1a0a8a3SJed Brown Mat Ad = PETSC_NULL, adj; 1129b1a0a8a3SJed Brown IS ispart,isnumb,*is; 1130b1a0a8a3SJed Brown PetscErrorCode ierr; 1131b1a0a8a3SJed Brown 1132b1a0a8a3SJed Brown PetscFunctionBegin; 1133b1a0a8a3SJed Brown PetscValidHeaderSpecific(A,MAT_CLASSID,1); 1134b1a0a8a3SJed Brown PetscValidPointer(outis,3); 1135b1a0a8a3SJed Brown if (n < 1) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of local blocks must be > 0, n = %D",n); 1136b1a0a8a3SJed Brown 1137b1a0a8a3SJed Brown /* Get prefix, row distribution, and block size */ 1138b1a0a8a3SJed Brown ierr = MatGetOptionsPrefix(A,&prefix);CHKERRQ(ierr); 1139b1a0a8a3SJed Brown ierr = MatGetOwnershipRange(A,&rstart,&rend);CHKERRQ(ierr); 1140b1a0a8a3SJed Brown ierr = MatGetBlockSize(A,&bs);CHKERRQ(ierr); 1141b1a0a8a3SJed 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); 1142b1a0a8a3SJed Brown 1143b1a0a8a3SJed Brown /* Get diagonal block from matrix if possible */ 1144b1a0a8a3SJed Brown ierr = MPI_Comm_size(((PetscObject)A)->comm,&size);CHKERRQ(ierr); 1145b1a0a8a3SJed Brown ierr = PetscObjectQueryFunction((PetscObject)A,"MatGetDiagonalBlock_C",(void (**)(void))&f);CHKERRQ(ierr); 1146b1a0a8a3SJed Brown if (f) { 1147b1a0a8a3SJed Brown ierr = (*f)(A,&iscopy,MAT_INITIAL_MATRIX,&Ad);CHKERRQ(ierr); 1148b1a0a8a3SJed Brown } else if (size == 1) { 1149b1a0a8a3SJed Brown iscopy = PETSC_FALSE; Ad = A; 1150b1a0a8a3SJed Brown } else { 1151b1a0a8a3SJed Brown iscopy = PETSC_FALSE; Ad = PETSC_NULL; 1152b1a0a8a3SJed Brown } 1153b1a0a8a3SJed Brown if (Ad) { 1154b1a0a8a3SJed Brown ierr = PetscTypeCompare((PetscObject)Ad,MATSEQBAIJ,&isbaij);CHKERRQ(ierr); 1155b1a0a8a3SJed Brown if (!isbaij) {ierr = PetscTypeCompare((PetscObject)Ad,MATSEQSBAIJ,&isbaij);CHKERRQ(ierr);} 1156b1a0a8a3SJed Brown } 1157b1a0a8a3SJed Brown if (Ad && n > 1) { 1158b1a0a8a3SJed Brown PetscBool match,done; 1159b1a0a8a3SJed Brown /* Try to setup a good matrix partitioning if available */ 1160b1a0a8a3SJed Brown ierr = MatPartitioningCreate(PETSC_COMM_SELF,&mpart);CHKERRQ(ierr); 1161b1a0a8a3SJed Brown ierr = PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix);CHKERRQ(ierr); 1162b1a0a8a3SJed Brown ierr = MatPartitioningSetFromOptions(mpart);CHKERRQ(ierr); 1163b1a0a8a3SJed Brown ierr = PetscTypeCompare((PetscObject)mpart,MATPARTITIONINGCURRENT,&match);CHKERRQ(ierr); 1164b1a0a8a3SJed Brown if (!match) { 1165b1a0a8a3SJed Brown ierr = PetscTypeCompare((PetscObject)mpart,MATPARTITIONINGSQUARE,&match);CHKERRQ(ierr); 1166b1a0a8a3SJed Brown } 1167b1a0a8a3SJed Brown if (!match) { /* assume a "good" partitioner is available */ 1168b1a0a8a3SJed Brown PetscInt na,*ia,*ja; 1169b1a0a8a3SJed Brown ierr = MatGetRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr); 1170b1a0a8a3SJed Brown if (done) { 1171b1a0a8a3SJed Brown /* Build adjacency matrix by hand. Unfortunately a call to 1172b1a0a8a3SJed Brown MatConvert(Ad,MATMPIADJ,MAT_INITIAL_MATRIX,&adj) will 1173b1a0a8a3SJed Brown remove the block-aij structure and we cannot expect 1174b1a0a8a3SJed Brown MatPartitioning to split vertices as we need */ 1175b1a0a8a3SJed Brown PetscInt i,j,*row,len,nnz,cnt,*iia=0,*jja=0; 1176b1a0a8a3SJed Brown nnz = 0; 1177b1a0a8a3SJed Brown for (i=0; i<na; i++) { /* count number of nonzeros */ 1178b1a0a8a3SJed Brown len = ia[i+1] - ia[i]; 1179b1a0a8a3SJed Brown row = ja + ia[i]; 1180b1a0a8a3SJed Brown for (j=0; j<len; j++) { 1181b1a0a8a3SJed Brown if (row[j] == i) { /* don't count diagonal */ 1182b1a0a8a3SJed Brown len--; break; 1183b1a0a8a3SJed Brown } 1184b1a0a8a3SJed Brown } 1185b1a0a8a3SJed Brown nnz += len; 1186b1a0a8a3SJed Brown } 1187b1a0a8a3SJed Brown ierr = PetscMalloc((na+1)*sizeof(PetscInt),&iia);CHKERRQ(ierr); 1188b1a0a8a3SJed Brown ierr = PetscMalloc((nnz)*sizeof(PetscInt),&jja);CHKERRQ(ierr); 1189b1a0a8a3SJed Brown nnz = 0; 1190b1a0a8a3SJed Brown iia[0] = 0; 1191b1a0a8a3SJed Brown for (i=0; i<na; i++) { /* fill adjacency */ 1192b1a0a8a3SJed Brown cnt = 0; 1193b1a0a8a3SJed Brown len = ia[i+1] - ia[i]; 1194b1a0a8a3SJed Brown row = ja + ia[i]; 1195b1a0a8a3SJed Brown for (j=0; j<len; j++) { 1196b1a0a8a3SJed Brown if (row[j] != i) { /* if not diagonal */ 1197b1a0a8a3SJed Brown jja[nnz+cnt++] = row[j]; 1198b1a0a8a3SJed Brown } 1199b1a0a8a3SJed Brown } 1200b1a0a8a3SJed Brown nnz += cnt; 1201b1a0a8a3SJed Brown iia[i+1] = nnz; 1202b1a0a8a3SJed Brown } 1203b1a0a8a3SJed Brown /* Partitioning of the adjacency matrix */ 1204b1a0a8a3SJed Brown ierr = MatCreateMPIAdj(PETSC_COMM_SELF,na,na,iia,jja,PETSC_NULL,&adj);CHKERRQ(ierr); 1205b1a0a8a3SJed Brown ierr = MatPartitioningSetAdjacency(mpart,adj);CHKERRQ(ierr); 1206b1a0a8a3SJed Brown ierr = MatPartitioningSetNParts(mpart,n);CHKERRQ(ierr); 1207b1a0a8a3SJed Brown ierr = MatPartitioningApply(mpart,&ispart);CHKERRQ(ierr); 1208b1a0a8a3SJed Brown ierr = ISPartitioningToNumbering(ispart,&isnumb);CHKERRQ(ierr); 1209b1a0a8a3SJed Brown ierr = MatDestroy(adj);CHKERRQ(ierr); 1210b1a0a8a3SJed Brown foundpart = PETSC_TRUE; 1211b1a0a8a3SJed Brown } 1212b1a0a8a3SJed Brown ierr = MatRestoreRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done);CHKERRQ(ierr); 1213b1a0a8a3SJed Brown } 1214b1a0a8a3SJed Brown ierr = MatPartitioningDestroy(mpart);CHKERRQ(ierr); 1215b1a0a8a3SJed Brown } 1216b1a0a8a3SJed Brown if (iscopy) {ierr = MatDestroy(Ad);CHKERRQ(ierr);} 1217b1a0a8a3SJed Brown 1218b1a0a8a3SJed Brown ierr = PetscMalloc(n*sizeof(IS),&is);CHKERRQ(ierr); 1219b1a0a8a3SJed Brown *outis = is; 1220b1a0a8a3SJed Brown 1221b1a0a8a3SJed Brown if (!foundpart) { 1222b1a0a8a3SJed Brown 1223b1a0a8a3SJed Brown /* Partitioning by contiguous chunks of rows */ 1224b1a0a8a3SJed Brown 1225b1a0a8a3SJed Brown PetscInt mbs = (rend-rstart)/bs; 1226b1a0a8a3SJed Brown PetscInt start = rstart; 1227b1a0a8a3SJed Brown for (i=0; i<n; i++) { 1228b1a0a8a3SJed Brown PetscInt count = (mbs/n + ((mbs % n) > i)) * bs; 1229b1a0a8a3SJed Brown ierr = ISCreateStride(PETSC_COMM_SELF,count,start,1,&is[i]);CHKERRQ(ierr); 1230b1a0a8a3SJed Brown start += count; 1231b1a0a8a3SJed Brown } 1232b1a0a8a3SJed Brown 1233b1a0a8a3SJed Brown } else { 1234b1a0a8a3SJed Brown 1235b1a0a8a3SJed Brown /* Partitioning by adjacency of diagonal block */ 1236b1a0a8a3SJed Brown 1237b1a0a8a3SJed Brown const PetscInt *numbering; 1238b1a0a8a3SJed Brown PetscInt *count,nidx,*indices,*newidx,start=0; 1239b1a0a8a3SJed Brown /* Get node count in each partition */ 1240b1a0a8a3SJed Brown ierr = PetscMalloc(n*sizeof(PetscInt),&count);CHKERRQ(ierr); 1241b1a0a8a3SJed Brown ierr = ISPartitioningCount(ispart,n,count);CHKERRQ(ierr); 1242b1a0a8a3SJed Brown if (isbaij && bs > 1) { /* adjust for the block-aij case */ 1243b1a0a8a3SJed Brown for (i=0; i<n; i++) count[i] *= bs; 1244b1a0a8a3SJed Brown } 1245b1a0a8a3SJed Brown /* Build indices from node numbering */ 1246b1a0a8a3SJed Brown ierr = ISGetLocalSize(isnumb,&nidx);CHKERRQ(ierr); 1247b1a0a8a3SJed Brown ierr = PetscMalloc(nidx*sizeof(PetscInt),&indices);CHKERRQ(ierr); 1248b1a0a8a3SJed Brown for (i=0; i<nidx; i++) indices[i] = i; /* needs to be initialized */ 1249b1a0a8a3SJed Brown ierr = ISGetIndices(isnumb,&numbering);CHKERRQ(ierr); 1250b1a0a8a3SJed Brown ierr = PetscSortIntWithPermutation(nidx,numbering,indices);CHKERRQ(ierr); 1251b1a0a8a3SJed Brown ierr = ISRestoreIndices(isnumb,&numbering);CHKERRQ(ierr); 1252b1a0a8a3SJed Brown if (isbaij && bs > 1) { /* adjust for the block-aij case */ 1253b1a0a8a3SJed Brown ierr = PetscMalloc(nidx*bs*sizeof(PetscInt),&newidx);CHKERRQ(ierr); 1254b1a0a8a3SJed Brown for (i=0; i<nidx; i++) 1255b1a0a8a3SJed Brown for (j=0; j<bs; j++) 1256b1a0a8a3SJed Brown newidx[i*bs+j] = indices[i]*bs + j; 1257b1a0a8a3SJed Brown ierr = PetscFree(indices);CHKERRQ(ierr); 1258b1a0a8a3SJed Brown nidx *= bs; 1259b1a0a8a3SJed Brown indices = newidx; 1260b1a0a8a3SJed Brown } 1261b1a0a8a3SJed Brown /* Shift to get global indices */ 1262b1a0a8a3SJed Brown for (i=0; i<nidx; i++) indices[i] += rstart; 1263b1a0a8a3SJed Brown 1264b1a0a8a3SJed Brown /* Build the index sets for each block */ 1265b1a0a8a3SJed Brown for (i=0; i<n; i++) { 1266b1a0a8a3SJed Brown ierr = ISCreateGeneral(PETSC_COMM_SELF,count[i],&indices[start],PETSC_COPY_VALUES,&is[i]);CHKERRQ(ierr); 1267b1a0a8a3SJed Brown ierr = ISSort(is[i]);CHKERRQ(ierr); 1268b1a0a8a3SJed Brown start += count[i]; 1269b1a0a8a3SJed Brown } 1270b1a0a8a3SJed Brown 1271b1a0a8a3SJed Brown ierr = PetscFree(count); 1272b1a0a8a3SJed Brown ierr = PetscFree(indices); 1273b1a0a8a3SJed Brown ierr = ISDestroy(isnumb);CHKERRQ(ierr); 1274b1a0a8a3SJed Brown ierr = ISDestroy(ispart);CHKERRQ(ierr); 1275b1a0a8a3SJed Brown } 1276b1a0a8a3SJed Brown 1277b1a0a8a3SJed Brown PetscFunctionReturn(0); 1278b1a0a8a3SJed Brown } 1279b1a0a8a3SJed Brown 1280b1a0a8a3SJed Brown #undef __FUNCT__ 1281f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMDestroySubdomains" 1282b1a0a8a3SJed Brown /*@C 1283f746d493SDmitry Karpeev PCGASMDestroySubdomains - Destroys the index sets created with 1284f746d493SDmitry Karpeev PCGASMCreateSubdomains(). Should be called after setting subdomains 1285f746d493SDmitry Karpeev with PCGASMSetLocalSubdomains(). 1286b1a0a8a3SJed Brown 1287b1a0a8a3SJed Brown Collective 1288b1a0a8a3SJed Brown 1289b1a0a8a3SJed Brown Input Parameters: 1290b1a0a8a3SJed Brown + n - the number of index sets 1291b1a0a8a3SJed Brown . is - the array of index sets 1292b1a0a8a3SJed Brown - is_local - the array of local index sets, can be PETSC_NULL 1293b1a0a8a3SJed Brown 1294b1a0a8a3SJed Brown Level: advanced 1295b1a0a8a3SJed Brown 1296f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, create, subdomains, unstructured grid 1297b1a0a8a3SJed Brown 1298f746d493SDmitry Karpeev .seealso: PCGASMCreateSubdomains(), PCGASMSetLocalSubdomains() 1299b1a0a8a3SJed Brown @*/ 13007087cfbeSBarry Smith PetscErrorCode PCGASMDestroySubdomains(PetscInt n, IS is[], IS is_local[]) 1301b1a0a8a3SJed Brown { 1302b1a0a8a3SJed Brown PetscInt i; 1303b1a0a8a3SJed Brown PetscErrorCode ierr; 1304b1a0a8a3SJed Brown PetscFunctionBegin; 1305b1a0a8a3SJed Brown if (n <= 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"n must be > 0: n = %D",n); 1306b1a0a8a3SJed Brown PetscValidPointer(is,2); 1307b1a0a8a3SJed Brown for (i=0; i<n; i++) { ierr = ISDestroy(is[i]);CHKERRQ(ierr); } 1308b1a0a8a3SJed Brown ierr = PetscFree(is);CHKERRQ(ierr); 1309b1a0a8a3SJed Brown if (is_local) { 1310b1a0a8a3SJed Brown PetscValidPointer(is_local,3); 1311b1a0a8a3SJed Brown for (i=0; i<n; i++) { ierr = ISDestroy(is_local[i]);CHKERRQ(ierr); } 1312b1a0a8a3SJed Brown ierr = PetscFree(is_local);CHKERRQ(ierr); 1313b1a0a8a3SJed Brown } 1314b1a0a8a3SJed Brown PetscFunctionReturn(0); 1315b1a0a8a3SJed Brown } 1316b1a0a8a3SJed Brown 1317af538404SDmitry Karpeev 1318af538404SDmitry Karpeev #define PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,xleft_loc,ylow_loc,xright_loc,yhigh_loc,n) \ 1319af538404SDmitry Karpeev { \ 1320af538404SDmitry Karpeev PetscInt first_row = first/M, last_row = last/M+1; \ 1321af538404SDmitry Karpeev /* \ 1322af538404SDmitry Karpeev Compute ylow_loc and yhigh_loc so that (ylow_loc,xleft) and (yhigh_loc,xright) are the corners \ 1323af538404SDmitry Karpeev of the bounding box of the intersection of the subdomain with the local ownership range (local \ 1324af538404SDmitry Karpeev subdomain). \ 1325af538404SDmitry Karpeev Also compute xleft_loc and xright_loc as the lower and upper bounds on the first and last rows \ 1326af538404SDmitry Karpeev of the intersection. \ 1327af538404SDmitry Karpeev */ \ 1328af538404SDmitry Karpeev /* ylow_loc is the grid row containing the first element of the local sumbdomain */ \ 1329eec7e3faSDmitry Karpeev *ylow_loc = PetscMax(first_row,ylow); \ 1330af538404SDmitry Karpeev /* xleft_loc is the offset of first element of the local subdomain within its grid row (might actually be outside the local subdomain) */ \ 1331eec7e3faSDmitry Karpeev *xleft_loc = *ylow_loc==first_row?PetscMax(first%M,xleft):xleft; \ 1332af538404SDmitry Karpeev /* yhigh_loc is the grid row above the last local subdomain element */ \ 1333eec7e3faSDmitry Karpeev *yhigh_loc = PetscMin(last_row,yhigh); \ 1334af538404SDmitry Karpeev /* xright is the offset of the end of the local subdomain within its grid row (might actually be outside the local subdomain) */ \ 1335eec7e3faSDmitry Karpeev *xright_loc = *yhigh_loc==last_row?PetscMin(xright,last%M):xright; \ 1336af538404SDmitry Karpeev /* Now compute the size of the local subdomain n. */ \ 1337c3518bceSDmitry Karpeev *n = 0; \ 1338eec7e3faSDmitry Karpeev if(*ylow_loc < *yhigh_loc) { \ 1339af538404SDmitry Karpeev PetscInt width = xright-xleft; \ 1340eec7e3faSDmitry Karpeev *n += width*(*yhigh_loc-*ylow_loc-1); \ 1341eec7e3faSDmitry Karpeev *n += PetscMin(PetscMax(*xright_loc-xleft,0),width); \ 1342eec7e3faSDmitry Karpeev *n -= PetscMin(PetscMax(*xleft_loc-xleft,0), width); \ 1343af538404SDmitry Karpeev }\ 1344af538404SDmitry Karpeev } 1345af538404SDmitry Karpeev 1346c3518bceSDmitry Karpeev 1347eec7e3faSDmitry Karpeev 1348eec7e3faSDmitry Karpeev #undef __FUNCT__ 1349f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMCreateSubdomains2D" 1350b1a0a8a3SJed Brown /*@ 1351f746d493SDmitry Karpeev PCGASMCreateSubdomains2D - Creates the index sets for the overlapping Schwarz 1352b1a0a8a3SJed Brown preconditioner for a two-dimensional problem on a regular grid. 1353b1a0a8a3SJed Brown 1354af538404SDmitry Karpeev Collective 1355b1a0a8a3SJed Brown 1356b1a0a8a3SJed Brown Input Parameters: 1357af538404SDmitry Karpeev + M, N - the global number of mesh points in the x and y directions 1358af538404SDmitry Karpeev . Mdomains, Ndomains - the global number of subdomains in the x and y directions 1359b1a0a8a3SJed Brown . dof - degrees of freedom per node 1360b1a0a8a3SJed Brown - overlap - overlap in mesh lines 1361b1a0a8a3SJed Brown 1362b1a0a8a3SJed Brown Output Parameters: 1363af538404SDmitry Karpeev + Nsub - the number of local subdomains created 1364b1a0a8a3SJed Brown . is - array of index sets defining overlapping (if overlap > 0) subdomains 1365b1a0a8a3SJed Brown - is_local - array of index sets defining non-overlapping subdomains 1366b1a0a8a3SJed Brown 1367b1a0a8a3SJed Brown Note: 1368b1a0a8a3SJed Brown Presently PCAMSCreateSubdomains2d() is valid only for sequential 1369b1a0a8a3SJed Brown preconditioners. More general related routines are 1370f746d493SDmitry Karpeev PCGASMSetTotalSubdomains() and PCGASMSetLocalSubdomains(). 1371b1a0a8a3SJed Brown 1372b1a0a8a3SJed Brown Level: advanced 1373b1a0a8a3SJed Brown 1374f746d493SDmitry Karpeev .keywords: PC, GASM, additive Schwarz, create, subdomains, 2D, regular grid 1375b1a0a8a3SJed Brown 1376f746d493SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetLocalSubdomains(), PCGASMGetSubKSP(), 1377f746d493SDmitry Karpeev PCGASMSetOverlap() 1378b1a0a8a3SJed Brown @*/ 13797087cfbeSBarry Smith PetscErrorCode PCGASMCreateSubdomains2D(PC pc, PetscInt M,PetscInt N,PetscInt Mdomains,PetscInt Ndomains,PetscInt dof,PetscInt overlap, PetscInt *nsub,IS **is,IS **is_local) 1380b1a0a8a3SJed Brown { 1381b1a0a8a3SJed Brown PetscErrorCode ierr; 13822bbb417fSDmitry Karpeev PetscMPIInt size, rank; 13832bbb417fSDmitry Karpeev PetscInt i, j; 13842bbb417fSDmitry Karpeev PetscInt maxheight, maxwidth; 13852bbb417fSDmitry Karpeev PetscInt xstart, xleft, xright, xleft_loc, xright_loc; 13862bbb417fSDmitry Karpeev PetscInt ystart, ylow, yhigh, ylow_loc, yhigh_loc; 1387eec7e3faSDmitry Karpeev PetscInt x[2][2], y[2][2], n[2]; 13882bbb417fSDmitry Karpeev PetscInt first, last; 13892bbb417fSDmitry Karpeev PetscInt nidx, *idx; 13902bbb417fSDmitry Karpeev PetscInt ii,jj,s,q,d; 1391af538404SDmitry Karpeev PetscInt k,kk; 13922bbb417fSDmitry Karpeev PetscMPIInt color; 13932bbb417fSDmitry Karpeev MPI_Comm comm, subcomm; 13945d38c402SBarry Smith IS **iis = 0; 1395d34fcf5fSBarry Smith 1396b1a0a8a3SJed Brown PetscFunctionBegin; 13972bbb417fSDmitry Karpeev ierr = PetscObjectGetComm((PetscObject)pc, &comm);CHKERRQ(ierr); 13982bbb417fSDmitry Karpeev ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 13992bbb417fSDmitry Karpeev ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 14002bbb417fSDmitry Karpeev ierr = MatGetOwnershipRange(pc->pmat, &first, &last);CHKERRQ(ierr); 1401d34fcf5fSBarry 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) " 14022bbb417fSDmitry Karpeev "does not respect the number of degrees of freedom per grid point %D", first, last, dof); 1403d34fcf5fSBarry Smith 1404af538404SDmitry Karpeev /* Determine the number of domains with nonzero intersections with the local ownership range. */ 14052bbb417fSDmitry Karpeev s = 0; 1406af538404SDmitry Karpeev ystart = 0; 1407af538404SDmitry Karpeev for (j=0; j<Ndomains; ++j) { 1408af538404SDmitry Karpeev maxheight = N/Ndomains + ((N % Ndomains) > j); /* Maximal height of subdomain */ 1409af538404SDmitry 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); 1410eec7e3faSDmitry Karpeev /* Vertical domain limits with an overlap. */ 1411eec7e3faSDmitry Karpeev ylow = PetscMax(ystart - overlap,0); 1412eec7e3faSDmitry Karpeev yhigh = PetscMin(ystart + maxheight + overlap,N); 1413b1a0a8a3SJed Brown xstart = 0; 1414af538404SDmitry Karpeev for (i=0; i<Mdomains; ++i) { 1415af538404SDmitry Karpeev maxwidth = M/Mdomains + ((M % Mdomains) > i); /* Maximal width of subdomain */ 1416af538404SDmitry 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); 1417eec7e3faSDmitry Karpeev /* Horizontal domain limits with an overlap. */ 1418eec7e3faSDmitry Karpeev xleft = PetscMax(xstart - overlap,0); 1419eec7e3faSDmitry Karpeev xright = PetscMin(xstart + maxwidth + overlap,M); 1420af538404SDmitry Karpeev /* 1421af538404SDmitry Karpeev Determine whether this subdomain intersects this processor's ownership range of pc->pmat. 1422af538404SDmitry Karpeev */ 1423c3518bceSDmitry Karpeev PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,(&xleft_loc),(&ylow_loc),(&xright_loc),(&yhigh_loc),(&nidx)); 1424af538404SDmitry Karpeev if(nidx) { 1425af538404SDmitry Karpeev ++s; 1426af538404SDmitry Karpeev } 1427af538404SDmitry Karpeev xstart += maxwidth; 1428af538404SDmitry Karpeev }/* for(i = 0; i < Mdomains; ++i) */ 1429af538404SDmitry Karpeev ystart += maxheight; 1430af538404SDmitry Karpeev }/* for(j = 0; j < Ndomains; ++j) */ 1431af538404SDmitry Karpeev /* Now we can allocate the necessary number of ISs. */ 1432af538404SDmitry Karpeev *nsub = s; 1433af538404SDmitry Karpeev ierr = PetscMalloc((*nsub)*sizeof(IS*),is);CHKERRQ(ierr); 1434af538404SDmitry Karpeev ierr = PetscMalloc((*nsub)*sizeof(IS*),is_local);CHKERRQ(ierr); 1435af538404SDmitry Karpeev s = 0; 1436af538404SDmitry Karpeev ystart = 0; 1437af538404SDmitry Karpeev for (j=0; j<Ndomains; ++j) { 1438af538404SDmitry Karpeev maxheight = N/Ndomains + ((N % Ndomains) > j); /* Maximal height of subdomain */ 1439af538404SDmitry 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); 1440af538404SDmitry Karpeev /* Vertical domain limits with an overlap. */ 1441eec7e3faSDmitry Karpeev y[0][0] = PetscMax(ystart - overlap,0); 1442eec7e3faSDmitry Karpeev y[0][1] = PetscMin(ystart + maxheight + overlap,N); 1443eec7e3faSDmitry Karpeev /* Vertical domain limits without an overlap. */ 1444eec7e3faSDmitry Karpeev y[1][0] = ystart; 1445eec7e3faSDmitry Karpeev y[1][1] = PetscMin(ystart + maxheight,N); 1446eec7e3faSDmitry Karpeev xstart = 0; 1447af538404SDmitry Karpeev for (i=0; i<Mdomains; ++i) { 1448af538404SDmitry Karpeev maxwidth = M/Mdomains + ((M % Mdomains) > i); /* Maximal width of subdomain */ 1449af538404SDmitry 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); 1450af538404SDmitry Karpeev /* Horizontal domain limits with an overlap. */ 1451eec7e3faSDmitry Karpeev x[0][0] = PetscMax(xstart - overlap,0); 1452eec7e3faSDmitry Karpeev x[0][1] = PetscMin(xstart + maxwidth + overlap,M); 1453eec7e3faSDmitry Karpeev /* Horizontal domain limits without an overlap. */ 1454eec7e3faSDmitry Karpeev x[1][0] = xstart; 1455eec7e3faSDmitry Karpeev x[1][1] = PetscMin(xstart+maxwidth,M); 14562bbb417fSDmitry Karpeev /* 14572bbb417fSDmitry Karpeev Determine whether this domain intersects this processor's ownership range of pc->pmat. 14582bbb417fSDmitry Karpeev Do this twice: first for the domains with overlaps, and once without. 14592bbb417fSDmitry Karpeev During the first pass create the subcommunicators, and use them on the second pass as well. 14602bbb417fSDmitry Karpeev */ 14612bbb417fSDmitry Karpeev for(q = 0; q < 2; ++q) { 14622bbb417fSDmitry Karpeev /* 14632bbb417fSDmitry Karpeev domain limits, (xleft, xright) and (ylow, yheigh) are adjusted 14642bbb417fSDmitry Karpeev according to whether the domain with an overlap or without is considered. 14652bbb417fSDmitry Karpeev */ 14662bbb417fSDmitry Karpeev xleft = x[q][0]; xright = x[q][1]; 14672bbb417fSDmitry Karpeev ylow = y[q][0]; yhigh = y[q][1]; 1468c3518bceSDmitry Karpeev PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,(&xleft_loc),(&ylow_loc),(&xright_loc),(&yhigh_loc),(&nidx)); 1469af538404SDmitry Karpeev nidx *= dof; 1470eec7e3faSDmitry Karpeev n[q] = nidx; 14712bbb417fSDmitry Karpeev /* 1472eec7e3faSDmitry Karpeev Based on the counted number of indices in the local domain *with an overlap*, 1473af538404SDmitry Karpeev construct a subcommunicator of all the processors supporting this domain. 1474eec7e3faSDmitry Karpeev Observe that a domain with an overlap might have nontrivial local support, 1475eec7e3faSDmitry Karpeev while the domain without an overlap might not. Hence, the decision to participate 1476eec7e3faSDmitry Karpeev in the subcommunicator must be based on the domain with an overlap. 14772bbb417fSDmitry Karpeev */ 14782bbb417fSDmitry Karpeev if (q == 0) { 14792bbb417fSDmitry Karpeev if(nidx) { 14802bbb417fSDmitry Karpeev color = 1; 1481d34fcf5fSBarry Smith } else { 14822bbb417fSDmitry Karpeev color = MPI_UNDEFINED; 14832bbb417fSDmitry Karpeev } 14842bbb417fSDmitry Karpeev ierr = MPI_Comm_split(comm, color, rank, &subcomm);CHKERRQ(ierr); 14852bbb417fSDmitry Karpeev } 1486af538404SDmitry Karpeev /* 1487eec7e3faSDmitry Karpeev Proceed only if the number of local indices *with an overlap* is nonzero. 1488af538404SDmitry Karpeev */ 1489eec7e3faSDmitry Karpeev if (n[0]) { 1490af538404SDmitry Karpeev if(q == 0) { 1491eec7e3faSDmitry Karpeev iis = is; 1492af538404SDmitry Karpeev } 1493af538404SDmitry Karpeev if (q == 1) { 1494af538404SDmitry Karpeev /* 1495eec7e3faSDmitry Karpeev The IS for the no-overlap subdomain shares a communicator with the overlapping domain. 1496af538404SDmitry Karpeev Moreover, if the overlap is zero, the two ISs are identical. 1497af538404SDmitry Karpeev */ 1498b1a0a8a3SJed Brown if (overlap == 0) { 1499eec7e3faSDmitry Karpeev (*is_local)[s] = (*is)[s]; 15002bbb417fSDmitry Karpeev ierr = PetscObjectReference((PetscObject)(*is)[s]);CHKERRQ(ierr); 15012bbb417fSDmitry Karpeev continue; 1502d34fcf5fSBarry Smith } else { 1503eec7e3faSDmitry Karpeev iis = is_local; 1504eec7e3faSDmitry Karpeev subcomm = ((PetscObject)(*is)[s])->comm; 15052bbb417fSDmitry Karpeev } 1506af538404SDmitry Karpeev }/* if(q == 1) */ 1507eec7e3faSDmitry Karpeev idx = PETSC_NULL; 15082bbb417fSDmitry Karpeev ierr = PetscMalloc(nidx*sizeof(PetscInt),&idx);CHKERRQ(ierr); 1509eec7e3faSDmitry Karpeev if(nidx) { 15102bbb417fSDmitry Karpeev k = 0; 15112bbb417fSDmitry Karpeev for (jj=ylow_loc; jj<yhigh_loc; ++jj) { 1512af538404SDmitry Karpeev PetscInt x0 = (jj==ylow_loc)?xleft_loc:xleft; 1513af538404SDmitry Karpeev PetscInt x1 = (jj==yhigh_loc-1)?xright_loc:xright; 1514af538404SDmitry Karpeev kk = dof*(M*jj + x0); 15152bbb417fSDmitry Karpeev for (ii=x0; ii<x1; ++ii) { 15162bbb417fSDmitry Karpeev for(d = 0; d < dof; ++d) { 15172bbb417fSDmitry Karpeev idx[k++] = kk++; 1518b1a0a8a3SJed Brown } 1519b1a0a8a3SJed Brown } 1520b1a0a8a3SJed Brown } 1521eec7e3faSDmitry Karpeev } 15222bbb417fSDmitry Karpeev ierr = ISCreateGeneral(subcomm,nidx,idx,PETSC_OWN_POINTER,(*iis)+s);CHKERRQ(ierr); 1523eec7e3faSDmitry Karpeev }/* if(n[0]) */ 15242bbb417fSDmitry Karpeev }/* for(q = 0; q < 2; ++q) */ 1525eec7e3faSDmitry Karpeev if(n[0]) { 15262bbb417fSDmitry Karpeev ++s; 1527b1a0a8a3SJed Brown } 1528af538404SDmitry Karpeev xstart += maxwidth; 1529af538404SDmitry Karpeev }/* for(i = 0; i < Mdomains; ++i) */ 15302bbb417fSDmitry Karpeev ystart += maxheight; 1531af538404SDmitry Karpeev }/* for(j = 0; j < Ndomains; ++j) */ 1532b1a0a8a3SJed Brown PetscFunctionReturn(0); 1533b1a0a8a3SJed Brown } 1534b1a0a8a3SJed Brown 1535b1a0a8a3SJed Brown #undef __FUNCT__ 1536f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMGetLocalSubdomains" 1537b1a0a8a3SJed Brown /*@C 1538f746d493SDmitry Karpeev PCGASMGetLocalSubdomains - Gets the local subdomains (for this processor 1539b1a0a8a3SJed Brown only) for the additive Schwarz preconditioner. 1540b1a0a8a3SJed Brown 1541b1a0a8a3SJed Brown Not Collective 1542b1a0a8a3SJed Brown 1543b1a0a8a3SJed Brown Input Parameter: 1544b1a0a8a3SJed Brown . pc - the preconditioner context 1545b1a0a8a3SJed Brown 1546b1a0a8a3SJed Brown Output Parameters: 1547b1a0a8a3SJed Brown + n - the number of subdomains for this processor (default value = 1) 1548b1a0a8a3SJed Brown . is - the index sets that define the subdomains for this processor 1549b1a0a8a3SJed Brown - is_local - the index sets that define the local part of the subdomains for this processor (can be PETSC_NULL) 1550b1a0a8a3SJed Brown 1551b1a0a8a3SJed Brown 1552b1a0a8a3SJed Brown Notes: 1553b1a0a8a3SJed Brown The IS numbering is in the parallel, global numbering of the vector. 1554b1a0a8a3SJed Brown 1555b1a0a8a3SJed Brown Level: advanced 1556b1a0a8a3SJed Brown 1557f746d493SDmitry Karpeev .keywords: PC, GASM, set, local, subdomains, additive Schwarz 1558b1a0a8a3SJed Brown 1559f746d493SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetOverlap(), PCGASMGetSubKSP(), 1560f746d493SDmitry Karpeev PCGASMCreateSubdomains2D(), PCGASMSetLocalSubdomains(), PCGASMGetLocalSubmatrices() 1561b1a0a8a3SJed Brown @*/ 15627087cfbeSBarry Smith PetscErrorCode PCGASMGetLocalSubdomains(PC pc,PetscInt *n,IS *is[],IS *is_local[]) 1563b1a0a8a3SJed Brown { 1564f746d493SDmitry Karpeev PC_GASM *osm; 1565b1a0a8a3SJed Brown PetscErrorCode ierr; 1566b1a0a8a3SJed Brown PetscBool match; 1567b1a0a8a3SJed Brown 1568b1a0a8a3SJed Brown PetscFunctionBegin; 1569b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1570b1a0a8a3SJed Brown PetscValidIntPointer(n,2); 1571b1a0a8a3SJed Brown if (is) PetscValidPointer(is,3); 1572f746d493SDmitry Karpeev ierr = PetscTypeCompare((PetscObject)pc,PCGASM,&match);CHKERRQ(ierr); 1573b1a0a8a3SJed Brown if (!match) { 1574b1a0a8a3SJed Brown if (n) *n = 0; 1575b1a0a8a3SJed Brown if (is) *is = PETSC_NULL; 1576b1a0a8a3SJed Brown } else { 1577f746d493SDmitry Karpeev osm = (PC_GASM*)pc->data; 1578ab3e923fSDmitry Karpeev if (n) *n = osm->n; 1579b1a0a8a3SJed Brown if (is) *is = osm->is; 1580b1a0a8a3SJed Brown if (is_local) *is_local = osm->is_local; 1581b1a0a8a3SJed Brown } 1582b1a0a8a3SJed Brown PetscFunctionReturn(0); 1583b1a0a8a3SJed Brown } 1584b1a0a8a3SJed Brown 1585b1a0a8a3SJed Brown #undef __FUNCT__ 1586f746d493SDmitry Karpeev #define __FUNCT__ "PCGASMGetLocalSubmatrices" 1587b1a0a8a3SJed Brown /*@C 1588f746d493SDmitry Karpeev PCGASMGetLocalSubmatrices - Gets the local submatrices (for this processor 1589b1a0a8a3SJed Brown only) for the additive Schwarz preconditioner. 1590b1a0a8a3SJed Brown 1591b1a0a8a3SJed Brown Not Collective 1592b1a0a8a3SJed Brown 1593b1a0a8a3SJed Brown Input Parameter: 1594b1a0a8a3SJed Brown . pc - the preconditioner context 1595b1a0a8a3SJed Brown 1596b1a0a8a3SJed Brown Output Parameters: 1597b1a0a8a3SJed Brown + n - the number of matrices for this processor (default value = 1) 1598b1a0a8a3SJed Brown - mat - the matrices 1599b1a0a8a3SJed Brown 1600b1a0a8a3SJed Brown 1601b1a0a8a3SJed Brown Level: advanced 1602b1a0a8a3SJed Brown 1603f746d493SDmitry Karpeev .keywords: PC, GASM, set, local, subdomains, additive Schwarz, block Jacobi 1604b1a0a8a3SJed Brown 1605f746d493SDmitry Karpeev .seealso: PCGASMSetTotalSubdomains(), PCGASMSetOverlap(), PCGASMGetSubKSP(), 1606f746d493SDmitry Karpeev PCGASMCreateSubdomains2D(), PCGASMSetLocalSubdomains(), PCGASMGetLocalSubdomains() 1607b1a0a8a3SJed Brown @*/ 16087087cfbeSBarry Smith PetscErrorCode PCGASMGetLocalSubmatrices(PC pc,PetscInt *n,Mat *mat[]) 1609b1a0a8a3SJed Brown { 1610f746d493SDmitry Karpeev PC_GASM *osm; 1611b1a0a8a3SJed Brown PetscErrorCode ierr; 1612b1a0a8a3SJed Brown PetscBool match; 1613b1a0a8a3SJed Brown 1614b1a0a8a3SJed Brown PetscFunctionBegin; 1615b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1616b1a0a8a3SJed Brown PetscValidIntPointer(n,2); 1617b1a0a8a3SJed Brown if (mat) PetscValidPointer(mat,3); 1618b1a0a8a3SJed Brown if (!pc->setupcalled) SETERRQ(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must call after KSPSetUP() or PCSetUp()."); 1619f746d493SDmitry Karpeev ierr = PetscTypeCompare((PetscObject)pc,PCGASM,&match);CHKERRQ(ierr); 1620b1a0a8a3SJed Brown if (!match) { 1621b1a0a8a3SJed Brown if (n) *n = 0; 1622b1a0a8a3SJed Brown if (mat) *mat = PETSC_NULL; 1623b1a0a8a3SJed Brown } else { 1624f746d493SDmitry Karpeev osm = (PC_GASM*)pc->data; 1625ab3e923fSDmitry Karpeev if (n) *n = osm->n; 1626b1a0a8a3SJed Brown if (mat) *mat = osm->pmat; 1627b1a0a8a3SJed Brown } 1628b1a0a8a3SJed Brown PetscFunctionReturn(0); 1629b1a0a8a3SJed Brown } 1630