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