1b1a0a8a3SJed Brown /* 2f746d493SDmitry Karpeev This file defines an "generalized" additive Schwarz preconditioner for any Mat implementation. 36a4f0f73SDmitry Karpeev In this version each processor may intersect multiple subdomains and any subdomain may 46a4f0f73SDmitry Karpeev intersect multiple processors. Intersections of subdomains with processors are called *local 56a4f0f73SDmitry Karpeev subdomains*. 6b1a0a8a3SJed Brown 78f3b4b4dSDmitry Karpeev N - total number of distinct global subdomains (set explicitly in PCGASMSetTotalSubdomains() or implicitly PCGASMSetSubdomains() and then calculated in PCSetUp_GASM()) 86a4f0f73SDmitry Karpeev n - actual number of local subdomains on this processor (set in PCGASMSetSubdomains() or calculated in PCGASMSetTotalSubdomains()) 98f3b4b4dSDmitry Karpeev nmax - maximum number of local subdomains per processor (calculated in PCSetUp_GASM()) 10b1a0a8a3SJed Brown */ 11af0996ceSBarry Smith #include <petsc/private/pcimpl.h> /*I "petscpc.h" I*/ 121e25c274SJed Brown #include <petscdm.h> 13b1a0a8a3SJed Brown 14b1a0a8a3SJed Brown typedef struct { 15f746d493SDmitry Karpeev PetscInt N,n,nmax; 16b1a0a8a3SJed Brown PetscInt overlap; /* overlap requested by user */ 178f3b4b4dSDmitry Karpeev PCGASMType type; /* use reduced interpolation, restriction or both */ 188f3b4b4dSDmitry Karpeev PetscBool type_set; /* if user set this value (so won't change it for symmetric problems) */ 198f3b4b4dSDmitry Karpeev PetscBool same_subdomain_solvers; /* flag indicating whether all local solvers are same */ 208f3b4b4dSDmitry Karpeev PetscBool sort_indices; /* flag to sort subdomain indices */ 218f3b4b4dSDmitry Karpeev PetscBool user_subdomains; /* whether the user set explicit subdomain index sets -- keep them on PCReset() */ 228f3b4b4dSDmitry Karpeev PetscBool dm_subdomains; /* whether DM is allowed to define subdomains */ 23ea91fabdSFande Kong PetscBool hierarchicalpartitioning; 248f3b4b4dSDmitry Karpeev IS *ois; /* index sets that define the outer (conceptually, overlapping) subdomains */ 258f3b4b4dSDmitry Karpeev IS *iis; /* index sets that define the inner (conceptually, nonoverlapping) subdomains */ 268f3b4b4dSDmitry Karpeev KSP *ksp; /* linear solvers for each subdomain */ 278f3b4b4dSDmitry Karpeev Mat *pmat; /* subdomain block matrices */ 28f746d493SDmitry Karpeev Vec gx,gy; /* Merged work vectors */ 29f746d493SDmitry Karpeev Vec *x,*y; /* Split work vectors; storage aliases pieces of storage of the above merged vectors. */ 306a4f0f73SDmitry Karpeev VecScatter gorestriction; /* merged restriction to disjoint union of outer subdomains */ 316a4f0f73SDmitry Karpeev VecScatter girestriction; /* merged restriction to disjoint union of inner subdomains */ 32ea91fabdSFande Kong VecScatter pctoouter; 33ea91fabdSFande Kong IS permutationIS; 34ea91fabdSFande Kong Mat permutationP; 35ea91fabdSFande Kong Mat pcmat; 36ea91fabdSFande Kong Vec pcx,pcy; 37f746d493SDmitry Karpeev } PC_GASM; 38b1a0a8a3SJed Brown 398f3b4b4dSDmitry Karpeev static PetscErrorCode PCGASMComputeGlobalSubdomainNumbering_Private(PC pc,PetscInt **numbering,PetscInt **permutation) 408f3b4b4dSDmitry Karpeev { 418f3b4b4dSDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 428f3b4b4dSDmitry Karpeev PetscInt i; 438f3b4b4dSDmitry Karpeev 448f3b4b4dSDmitry Karpeev PetscFunctionBegin; 458f3b4b4dSDmitry Karpeev /* Determine the number of globally-distinct subdomains and compute a global numbering for them. */ 469566063dSJacob Faibussowitsch PetscCall(PetscMalloc2(osm->n,numbering,osm->n,permutation)); 479566063dSJacob Faibussowitsch PetscCall(PetscObjectsListGetGlobalNumbering(PetscObjectComm((PetscObject)pc),osm->n,(PetscObject*)osm->iis,NULL,*numbering)); 488f3b4b4dSDmitry Karpeev for (i = 0; i < osm->n; ++i) (*permutation)[i] = i; 499566063dSJacob Faibussowitsch PetscCall(PetscSortIntWithPermutation(osm->n,*numbering,*permutation)); 508f3b4b4dSDmitry Karpeev PetscFunctionReturn(0); 518f3b4b4dSDmitry Karpeev } 528f3b4b4dSDmitry Karpeev 5306b43e7eSDmitry Karpeev static PetscErrorCode PCGASMSubdomainView_Private(PC pc, PetscInt i, PetscViewer viewer) 54af538404SDmitry Karpeev { 55af538404SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 5606b43e7eSDmitry Karpeev PetscInt j,nidx; 57af538404SDmitry Karpeev const PetscInt *idx; 5806b43e7eSDmitry Karpeev PetscViewer sviewer; 59af538404SDmitry Karpeev char *cidx; 60af538404SDmitry Karpeev 61af538404SDmitry Karpeev PetscFunctionBegin; 622472a847SBarry Smith PetscCheck(i >= 0 && i < osm->n,PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONG, "Invalid subdomain %" PetscInt_FMT ": must nonnegative and less than %" PetscInt_FMT, i, osm->n); 634bde246dSDmitry Karpeev /* Inner subdomains. */ 649566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(osm->iis[i], &nidx)); 654bde246dSDmitry Karpeev /* 664bde246dSDmitry Karpeev No more than 15 characters per index plus a space. 674bde246dSDmitry Karpeev PetscViewerStringSPrintf requires a string of size at least 2, so use (nidx+1) instead of nidx, 684bde246dSDmitry Karpeev in case nidx == 0. That will take care of the space for the trailing '\0' as well. 694bde246dSDmitry Karpeev For nidx == 0, the whole string 16 '\0'. 704bde246dSDmitry Karpeev */ 7136a9e3b9SBarry Smith #define len 16*(nidx+1)+1 729566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(len, &cidx)); 739566063dSJacob Faibussowitsch PetscCall(PetscViewerStringOpen(PETSC_COMM_SELF, cidx, len, &sviewer)); 7436a9e3b9SBarry Smith #undef len 759566063dSJacob Faibussowitsch PetscCall(ISGetIndices(osm->iis[i], &idx)); 764bde246dSDmitry Karpeev for (j = 0; j < nidx; ++j) { 7763a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerStringSPrintf(sviewer, "%" PetscInt_FMT " ", idx[j])); 784bde246dSDmitry Karpeev } 799566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(osm->iis[i],&idx)); 809566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&sviewer)); 819566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Inner subdomain:\n")); 829566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 839566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 849566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "%s", cidx)); 859566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 869566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 879566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "\n")); 889566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 899566063dSJacob Faibussowitsch PetscCall(PetscFree(cidx)); 904bde246dSDmitry Karpeev /* Outer subdomains. */ 919566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(osm->ois[i], &nidx)); 92eec7e3faSDmitry Karpeev /* 93eec7e3faSDmitry Karpeev No more than 15 characters per index plus a space. 94eec7e3faSDmitry Karpeev PetscViewerStringSPrintf requires a string of size at least 2, so use (nidx+1) instead of nidx, 95eec7e3faSDmitry Karpeev in case nidx == 0. That will take care of the space for the trailing '\0' as well. 96eec7e3faSDmitry Karpeev For nidx == 0, the whole string 16 '\0'. 97eec7e3faSDmitry Karpeev */ 9836a9e3b9SBarry Smith #define len 16*(nidx+1)+1 999566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(len, &cidx)); 1009566063dSJacob Faibussowitsch PetscCall(PetscViewerStringOpen(PETSC_COMM_SELF, cidx, len, &sviewer)); 10136a9e3b9SBarry Smith #undef len 1029566063dSJacob Faibussowitsch PetscCall(ISGetIndices(osm->ois[i], &idx)); 103af538404SDmitry Karpeev for (j = 0; j < nidx; ++j) { 10463a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerStringSPrintf(sviewer,"%" PetscInt_FMT " ", idx[j])); 105af538404SDmitry Karpeev } 1069566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&sviewer)); 1079566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(osm->ois[i],&idx)); 1089566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Outer subdomain:\n")); 1099566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 1109566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 1119566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "%s", cidx)); 1129566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 1139566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 1149566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "\n")); 1159566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 1169566063dSJacob Faibussowitsch PetscCall(PetscFree(cidx)); 11706b43e7eSDmitry Karpeev PetscFunctionReturn(0); 11806b43e7eSDmitry Karpeev } 11906b43e7eSDmitry Karpeev 12006b43e7eSDmitry Karpeev static PetscErrorCode PCGASMPrintSubdomains(PC pc) 12106b43e7eSDmitry Karpeev { 12206b43e7eSDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 12306b43e7eSDmitry Karpeev const char *prefix; 12406b43e7eSDmitry Karpeev char fname[PETSC_MAX_PATH_LEN+1]; 1258f3b4b4dSDmitry Karpeev PetscInt l, d, count; 1269140fbc5SPierre Jolivet PetscBool found; 1270298fd71SBarry Smith PetscViewer viewer, sviewer = NULL; 1288f3b4b4dSDmitry Karpeev PetscInt *numbering,*permutation;/* global numbering of locally-supported subdomains and the permutation from the local ordering */ 12906b43e7eSDmitry Karpeev 13006b43e7eSDmitry Karpeev PetscFunctionBegin; 1319566063dSJacob Faibussowitsch PetscCall(PCGetOptionsPrefix(pc,&prefix)); 1329566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL,prefix,"-pc_gasm_print_subdomains",&found)); 1339140fbc5SPierre Jolivet if (!found) PetscFunctionReturn(0); 1349566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetString(NULL,prefix,"-pc_gasm_print_subdomains",fname,sizeof(fname),&found)); 1359566063dSJacob Faibussowitsch if (!found) PetscCall(PetscStrcpy(fname,"stdout")); 1369566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIOpen(PetscObjectComm((PetscObject)pc),fname,&viewer)); 1374bde246dSDmitry Karpeev /* 1384bde246dSDmitry Karpeev Make sure the viewer has a name. Otherwise this may cause a deadlock or other weird errors when creating a subcomm viewer: 1394bde246dSDmitry Karpeev the subcomm viewer will attempt to inherit the viewer's name, which, if not set, will be constructed collectively on the comm. 1404bde246dSDmitry Karpeev */ 1419566063dSJacob Faibussowitsch PetscCall(PetscObjectName((PetscObject)viewer)); 1424bde246dSDmitry Karpeev l = 0; 1439566063dSJacob Faibussowitsch PetscCall(PCGASMComputeGlobalSubdomainNumbering_Private(pc,&numbering,&permutation)); 1448f3b4b4dSDmitry Karpeev for (count = 0; count < osm->N; ++count) { 1454bde246dSDmitry Karpeev /* Now let subdomains go one at a time in the global numbering order and print their subdomain/solver info. */ 1464bde246dSDmitry Karpeev if (l<osm->n) { 1474bde246dSDmitry Karpeev d = permutation[l]; /* d is the local number of the l-th smallest (in the global ordering) among the locally supported subdomains */ 1484bde246dSDmitry Karpeev if (numbering[d] == count) { 1499566063dSJacob Faibussowitsch PetscCall(PetscViewerGetSubViewer(viewer,((PetscObject)osm->ois[d])->comm, &sviewer)); 1509566063dSJacob Faibussowitsch PetscCall(PCGASMSubdomainView_Private(pc,d,sviewer)); 1519566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(viewer,((PetscObject)osm->ois[d])->comm, &sviewer)); 1524bde246dSDmitry Karpeev ++l; 153af538404SDmitry Karpeev } 1544bde246dSDmitry Karpeev } 1559566063dSJacob Faibussowitsch PetscCallMPI(MPI_Barrier(PetscObjectComm((PetscObject)pc))); 1564bde246dSDmitry Karpeev } 1579566063dSJacob Faibussowitsch PetscCall(PetscFree2(numbering,permutation)); 1589566063dSJacob Faibussowitsch PetscCall(PetscViewerDestroy(&viewer)); 159af538404SDmitry Karpeev PetscFunctionReturn(0); 160af538404SDmitry Karpeev } 161af538404SDmitry Karpeev 162f746d493SDmitry Karpeev static PetscErrorCode PCView_GASM(PC pc,PetscViewer viewer) 163b1a0a8a3SJed Brown { 164f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 165af538404SDmitry Karpeev const char *prefix; 166af538404SDmitry Karpeev PetscMPIInt rank, size; 1678f3b4b4dSDmitry Karpeev PetscInt bsz; 16806b43e7eSDmitry Karpeev PetscBool iascii,view_subdomains=PETSC_FALSE; 169b1a0a8a3SJed Brown PetscViewer sviewer; 1708f3b4b4dSDmitry Karpeev PetscInt count, l; 1716a4f0f73SDmitry Karpeev char overlap[256] = "user-defined overlap"; 1726a4f0f73SDmitry Karpeev char gsubdomains[256] = "unknown total number of subdomains"; 17306b43e7eSDmitry Karpeev char msubdomains[256] = "unknown max number of local subdomains"; 1748f3b4b4dSDmitry Karpeev PetscInt *numbering,*permutation;/* global numbering of locally-supported subdomains and the permutation from the local ordering */ 17511aeaf0aSBarry Smith 176b1a0a8a3SJed Brown PetscFunctionBegin; 1779566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size)); 1789566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank)); 17906b43e7eSDmitry Karpeev 18006b43e7eSDmitry Karpeev if (osm->overlap >= 0) { 18163a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(overlap,sizeof(overlap),"requested amount of overlap = %" PetscInt_FMT,osm->overlap)); 18206b43e7eSDmitry Karpeev } 1838f3b4b4dSDmitry Karpeev if (osm->N != PETSC_DETERMINE) { 18463a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(gsubdomains, sizeof(gsubdomains), "total number of subdomains = %" PetscInt_FMT,osm->N)); 18506b43e7eSDmitry Karpeev } 1868f3b4b4dSDmitry Karpeev if (osm->nmax != PETSC_DETERMINE) { 18763a3b9bcSJacob Faibussowitsch PetscCall(PetscSNPrintf(msubdomains,sizeof(msubdomains),"max number of local subdomains = %" PetscInt_FMT,osm->nmax)); 18806b43e7eSDmitry Karpeev } 18906b43e7eSDmitry Karpeev 1909566063dSJacob Faibussowitsch PetscCall(PCGetOptionsPrefix(pc,&prefix)); 1919566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL,prefix,"-pc_gasm_view_subdomains",&view_subdomains,NULL)); 19206b43e7eSDmitry Karpeev 1939566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii)); 194b1a0a8a3SJed Brown if (iascii) { 19506b43e7eSDmitry Karpeev /* 19606b43e7eSDmitry Karpeev Make sure the viewer has a name. Otherwise this may cause a deadlock when creating a subcomm viewer: 19706b43e7eSDmitry Karpeev the subcomm viewer will attempt to inherit the viewer's name, which, if not set, will be constructed 19806b43e7eSDmitry Karpeev collectively on the comm. 19906b43e7eSDmitry Karpeev */ 2009566063dSJacob Faibussowitsch PetscCall(PetscObjectName((PetscObject)viewer)); 2019566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Restriction/interpolation type: %s\n",PCGASMTypes[osm->type])); 2029566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," %s\n",overlap)); 2039566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," %s\n",gsubdomains)); 2049566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," %s\n",msubdomains)); 2059566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushSynchronized(viewer)); 20663a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(viewer," [%d|%d] number of locally-supported subdomains = %" PetscInt_FMT "\n",rank,size,osm->n)); 2079566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 2089566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 2096a4f0f73SDmitry Karpeev /* Cannot take advantage of osm->same_subdomain_solvers without a global numbering of subdomains. */ 2109566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Subdomain solver info is as follows:\n")); 2119566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 2129566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," - - - - - - - - - - - - - - - - - -\n")); 21306b43e7eSDmitry Karpeev /* Make sure that everybody waits for the banner to be printed. */ 2149566063dSJacob Faibussowitsch PetscCallMPI(MPI_Barrier(PetscObjectComm((PetscObject)viewer))); 2154bde246dSDmitry Karpeev /* Now let subdomains go one at a time in the global numbering order and print their subdomain/solver info. */ 2169566063dSJacob Faibussowitsch PetscCall(PCGASMComputeGlobalSubdomainNumbering_Private(pc,&numbering,&permutation)); 21706b43e7eSDmitry Karpeev l = 0; 2188f3b4b4dSDmitry Karpeev for (count = 0; count < osm->N; ++count) { 21906b43e7eSDmitry Karpeev PetscMPIInt srank, ssize; 22006b43e7eSDmitry Karpeev if (l<osm->n) { 22106b43e7eSDmitry Karpeev PetscInt d = permutation[l]; /* d is the local number of the l-th smallest (in the global ordering) among the locally supported subdomains */ 22206b43e7eSDmitry Karpeev if (numbering[d] == count) { 2239566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(((PetscObject)osm->ois[d])->comm, &ssize)); 2249566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(((PetscObject)osm->ois[d])->comm, &srank)); 2259566063dSJacob Faibussowitsch PetscCall(PetscViewerGetSubViewer(viewer,((PetscObject)osm->ois[d])->comm, &sviewer)); 2269566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(osm->ois[d],&bsz)); 22763a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIISynchronizedPrintf(sviewer," [%d|%d] (subcomm [%d|%d]) local subdomain number %" PetscInt_FMT ", local size = %" PetscInt_FMT "\n",rank,size,srank,ssize,d,bsz)); 2289566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(sviewer)); 229*1baa6e33SBarry Smith if (view_subdomains) PetscCall(PCGASMSubdomainView_Private(pc,d,sviewer)); 2306a4f0f73SDmitry Karpeev if (!pc->setupcalled) { 2319566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(sviewer, " Solver not set up yet: PCSetUp() not yet called\n")); 2322fa5cd67SKarl Rupp } else { 2339566063dSJacob Faibussowitsch PetscCall(KSPView(osm->ksp[d],sviewer)); 2346a4f0f73SDmitry Karpeev } 2359566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(sviewer," - - - - - - - - - - - - - - - - - -\n")); 2369566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(sviewer)); 2379566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(viewer,((PetscObject)osm->ois[d])->comm, &sviewer)); 23806b43e7eSDmitry Karpeev ++l; 239b1a0a8a3SJed Brown } 24006b43e7eSDmitry Karpeev } 2419566063dSJacob Faibussowitsch PetscCallMPI(MPI_Barrier(PetscObjectComm((PetscObject)pc))); 242b1a0a8a3SJed Brown } 2439566063dSJacob Faibussowitsch PetscCall(PetscFree2(numbering,permutation)); 2449566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 2459566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 2469530cbd7SBarry Smith /* this line is needed to match the extra PetscViewerASCIIPushSynchronized() in PetscViewerGetSubViewer() */ 2479566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 2489530cbd7SBarry Smith } 249b1a0a8a3SJed Brown PetscFunctionReturn(0); 250b1a0a8a3SJed Brown } 251b1a0a8a3SJed Brown 2528f3b4b4dSDmitry Karpeev PETSC_INTERN PetscErrorCode PCGASMCreateLocalSubdomains(Mat A, PetscInt nloc, IS *iis[]); 253af538404SDmitry Karpeev 254ea91fabdSFande Kong PetscErrorCode PCGASMSetHierarchicalPartitioning(PC pc) 255ea91fabdSFande Kong { 256ea91fabdSFande Kong PC_GASM *osm = (PC_GASM*)pc->data; 257ea91fabdSFande Kong MatPartitioning part; 258ea91fabdSFande Kong MPI_Comm comm; 259ea91fabdSFande Kong PetscMPIInt size; 260ea91fabdSFande Kong PetscInt nlocalsubdomains,fromrows_localsize; 261ea91fabdSFande Kong IS partitioning,fromrows,isn; 262ea91fabdSFande Kong Vec outervec; 263ea91fabdSFande Kong 264ea91fabdSFande Kong PetscFunctionBegin; 2659566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)pc,&comm)); 2669566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm,&size)); 267ea91fabdSFande Kong /* we do not need a hierarchical partitioning when 268ea91fabdSFande Kong * the total number of subdomains is consistent with 269ea91fabdSFande Kong * the number of MPI tasks. 270ea91fabdSFande Kong * For the following cases, we do not need to use HP 271ea91fabdSFande Kong * */ 272670417b2SFande Kong if (osm->N==PETSC_DETERMINE || osm->N>=size || osm->N==1) PetscFunctionReturn(0); 27363a3b9bcSJacob Faibussowitsch PetscCheck(size%osm->N == 0,PETSC_COMM_WORLD,PETSC_ERR_ARG_INCOMP,"have to specify the total number of subdomains %" PetscInt_FMT " to be a factor of the number of processors %d ",osm->N,size); 274ea91fabdSFande Kong nlocalsubdomains = size/osm->N; 275ea91fabdSFande Kong osm->n = 1; 2769566063dSJacob Faibussowitsch PetscCall(MatPartitioningCreate(comm,&part)); 2779566063dSJacob Faibussowitsch PetscCall(MatPartitioningSetAdjacency(part,pc->pmat)); 2789566063dSJacob Faibussowitsch PetscCall(MatPartitioningSetType(part,MATPARTITIONINGHIERARCH)); 2799566063dSJacob Faibussowitsch PetscCall(MatPartitioningHierarchicalSetNcoarseparts(part,osm->N)); 2809566063dSJacob Faibussowitsch PetscCall(MatPartitioningHierarchicalSetNfineparts(part,nlocalsubdomains)); 2819566063dSJacob Faibussowitsch PetscCall(MatPartitioningSetFromOptions(part)); 282ea91fabdSFande Kong /* get new processor owner number of each vertex */ 2839566063dSJacob Faibussowitsch PetscCall(MatPartitioningApply(part,&partitioning)); 2849566063dSJacob Faibussowitsch PetscCall(ISBuildTwoSided(partitioning,NULL,&fromrows)); 2859566063dSJacob Faibussowitsch PetscCall(ISPartitioningToNumbering(partitioning,&isn)); 2869566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isn)); 2879566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(fromrows,&fromrows_localsize)); 2889566063dSJacob Faibussowitsch PetscCall(MatPartitioningDestroy(&part)); 2899566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(pc->pmat,&outervec,NULL)); 2909566063dSJacob Faibussowitsch PetscCall(VecCreateMPI(comm,fromrows_localsize,PETSC_DETERMINE,&(osm->pcx))); 2919566063dSJacob Faibussowitsch PetscCall(VecDuplicate(osm->pcx,&(osm->pcy))); 2929566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(osm->pcx,NULL,outervec,fromrows,&(osm->pctoouter))); 2939566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(pc->pmat,fromrows,fromrows,MAT_INITIAL_MATRIX,&(osm->permutationP))); 2949566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)fromrows)); 295ea91fabdSFande Kong osm->permutationIS = fromrows; 296ea91fabdSFande Kong osm->pcmat = pc->pmat; 2979566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)osm->permutationP)); 298ea91fabdSFande Kong pc->pmat = osm->permutationP; 2999566063dSJacob Faibussowitsch PetscCall(VecDestroy(&outervec)); 3009566063dSJacob Faibussowitsch PetscCall(ISDestroy(&fromrows)); 3019566063dSJacob Faibussowitsch PetscCall(ISDestroy(&partitioning)); 302ea91fabdSFande Kong osm->n = PETSC_DETERMINE; 303ea91fabdSFande Kong PetscFunctionReturn(0); 304ea91fabdSFande Kong } 305ea91fabdSFande Kong 306f746d493SDmitry Karpeev static PetscErrorCode PCSetUp_GASM(PC pc) 307b1a0a8a3SJed Brown { 308f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 309ea91fabdSFande Kong PetscInt i,nInnerIndices,nTotalInnerIndices; 310c10bc1a1SDmitry Karpeev PetscMPIInt rank, size; 311b1a0a8a3SJed Brown MatReuse scall = MAT_REUSE_MATRIX; 312b1a0a8a3SJed Brown KSP ksp; 313b1a0a8a3SJed Brown PC subpc; 314b1a0a8a3SJed Brown const char *prefix,*pprefix; 315f746d493SDmitry Karpeev Vec x,y; 3166a4f0f73SDmitry Karpeev PetscInt oni; /* Number of indices in the i-th local outer subdomain. */ 3176a4f0f73SDmitry Karpeev const PetscInt *oidxi; /* Indices from the i-th subdomain local outer subdomain. */ 3186a4f0f73SDmitry Karpeev PetscInt on; /* Number of indices in the disjoint union of local outer subdomains. */ 3196a4f0f73SDmitry Karpeev PetscInt *oidx; /* Indices in the disjoint union of local outer subdomains. */ 3206a4f0f73SDmitry Karpeev IS gois; /* Disjoint union the global indices of outer subdomains. */ 3216a4f0f73SDmitry Karpeev IS goid; /* Identity IS of the size of the disjoint union of outer subdomains. */ 322f746d493SDmitry Karpeev PetscScalar *gxarray, *gyarray; 323930d09c1SFande Kong PetscInt gostart; /* Start of locally-owned indices in the vectors -- osm->gx,osm->gy -- over the disjoint union of outer subdomains. */ 3248f3b4b4dSDmitry Karpeev PetscInt num_subdomains = 0; 3250298fd71SBarry Smith DM *subdomain_dm = NULL; 3268f3b4b4dSDmitry Karpeev char **subdomain_names = NULL; 327f771a274SFande Kong PetscInt *numbering; 3288f3b4b4dSDmitry Karpeev 329b1a0a8a3SJed Brown PetscFunctionBegin; 3309566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size)); 3319566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank)); 332b1a0a8a3SJed Brown if (!pc->setupcalled) { 333ea91fabdSFande Kong /* use a hierarchical partitioning */ 334*1baa6e33SBarry Smith if (osm->hierarchicalpartitioning) PetscCall(PCGASMSetHierarchicalPartitioning(pc)); 3358f3b4b4dSDmitry Karpeev if (osm->n == PETSC_DETERMINE) { 3368f3b4b4dSDmitry Karpeev if (osm->N != PETSC_DETERMINE) { 3378f3b4b4dSDmitry Karpeev /* No local subdomains given, but the desired number of total subdomains is known, so construct them accordingly. */ 3389566063dSJacob Faibussowitsch PetscCall(PCGASMCreateSubdomains(pc->pmat,osm->N,&osm->n,&osm->iis)); 3398f3b4b4dSDmitry Karpeev } else if (osm->dm_subdomains && pc->dm) { 3408f3b4b4dSDmitry Karpeev /* try pc->dm next, if allowed */ 3418f3b4b4dSDmitry Karpeev PetscInt d; 342a35b7b57SDmitry Karpeev IS *inner_subdomain_is, *outer_subdomain_is; 3439566063dSJacob Faibussowitsch PetscCall(DMCreateDomainDecomposition(pc->dm, &num_subdomains, &subdomain_names, &inner_subdomain_is, &outer_subdomain_is, &subdomain_dm)); 344*1baa6e33SBarry Smith if (num_subdomains) PetscCall(PCGASMSetSubdomains(pc, num_subdomains, inner_subdomain_is, outer_subdomain_is)); 345a35b7b57SDmitry Karpeev for (d = 0; d < num_subdomains; ++d) { 3469566063dSJacob Faibussowitsch if (inner_subdomain_is) PetscCall(ISDestroy(&inner_subdomain_is[d])); 3479566063dSJacob Faibussowitsch if (outer_subdomain_is) PetscCall(ISDestroy(&outer_subdomain_is[d])); 348fdc48646SDmitry Karpeev } 3499566063dSJacob Faibussowitsch PetscCall(PetscFree(inner_subdomain_is)); 3509566063dSJacob Faibussowitsch PetscCall(PetscFree(outer_subdomain_is)); 3518f3b4b4dSDmitry Karpeev } else { 3528f3b4b4dSDmitry Karpeev /* still no subdomains; use one per processor */ 35344fe18e5SDmitry Karpeev osm->nmax = osm->n = 1; 3549566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size)); 355f746d493SDmitry Karpeev osm->N = size; 3569566063dSJacob Faibussowitsch PetscCall(PCGASMCreateLocalSubdomains(pc->pmat,osm->n,&osm->iis)); 357fdc48646SDmitry Karpeev } 35806b43e7eSDmitry Karpeev } 359a35b7b57SDmitry Karpeev if (!osm->iis) { 360a35b7b57SDmitry Karpeev /* 3618f3b4b4dSDmitry Karpeev osm->n was set in PCGASMSetSubdomains(), but the actual subdomains have not been supplied. 3628f3b4b4dSDmitry Karpeev We create the requisite number of local inner subdomains and then expand them into 3638f3b4b4dSDmitry Karpeev out subdomains, if necessary. 364a35b7b57SDmitry Karpeev */ 3659566063dSJacob Faibussowitsch PetscCall(PCGASMCreateLocalSubdomains(pc->pmat,osm->n,&osm->iis)); 366a35b7b57SDmitry Karpeev } 3678f3b4b4dSDmitry Karpeev if (!osm->ois) { 3688f3b4b4dSDmitry Karpeev /* 3698f3b4b4dSDmitry Karpeev Initially make outer subdomains the same as inner subdomains. If nonzero additional overlap 3708f3b4b4dSDmitry Karpeev has been requested, copy the inner subdomains over so they can be modified. 3718f3b4b4dSDmitry Karpeev */ 3729566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(osm->n,&osm->ois)); 3738f3b4b4dSDmitry Karpeev for (i=0; i<osm->n; ++i) { 374ea91fabdSFande Kong if (osm->overlap > 0 && osm->N>1) { /* With positive overlap, osm->iis[i] will be modified */ 3759566063dSJacob Faibussowitsch PetscCall(ISDuplicate(osm->iis[i],(osm->ois)+i)); 3769566063dSJacob Faibussowitsch PetscCall(ISCopy(osm->iis[i],osm->ois[i])); 3778f3b4b4dSDmitry Karpeev } else { 3789566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)((osm->iis)[i]))); 3798f3b4b4dSDmitry Karpeev osm->ois[i] = osm->iis[i]; 3808f3b4b4dSDmitry Karpeev } 3818f3b4b4dSDmitry Karpeev } 382ea91fabdSFande Kong if (osm->overlap>0 && osm->N>1) { 3838f3b4b4dSDmitry Karpeev /* Extend the "overlapping" regions by a number of steps */ 3849566063dSJacob Faibussowitsch PetscCall(MatIncreaseOverlapSplit(pc->pmat,osm->n,osm->ois,osm->overlap)); 3858f3b4b4dSDmitry Karpeev } 386b1a0a8a3SJed Brown } 3876a4f0f73SDmitry Karpeev 3888f3b4b4dSDmitry Karpeev /* Now the subdomains are defined. Determine their global and max local numbers, if necessary. */ 3898f3b4b4dSDmitry Karpeev if (osm->nmax == PETSC_DETERMINE) { 3908f3b4b4dSDmitry Karpeev PetscMPIInt inwork,outwork; 3918f3b4b4dSDmitry Karpeev /* determine global number of subdomains and the max number of local subdomains */ 3928f3b4b4dSDmitry Karpeev inwork = osm->n; 3931c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(&inwork,&outwork,1,MPI_INT,MPI_MAX,PetscObjectComm((PetscObject)pc))); 3948f3b4b4dSDmitry Karpeev osm->nmax = outwork; 3958f3b4b4dSDmitry Karpeev } 3968f3b4b4dSDmitry Karpeev if (osm->N == PETSC_DETERMINE) { 3978f3b4b4dSDmitry Karpeev /* Determine the number of globally-distinct subdomains and compute a global numbering for them. */ 3989566063dSJacob Faibussowitsch PetscCall(PetscObjectsListGetGlobalNumbering(PetscObjectComm((PetscObject)pc),osm->n,(PetscObject*)osm->ois,&osm->N,NULL)); 3998f3b4b4dSDmitry Karpeev } 4008f3b4b4dSDmitry Karpeev 401b1a0a8a3SJed Brown if (osm->sort_indices) { 402f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 4039566063dSJacob Faibussowitsch PetscCall(ISSort(osm->ois[i])); 4049566063dSJacob Faibussowitsch PetscCall(ISSort(osm->iis[i])); 405b1a0a8a3SJed Brown } 406b1a0a8a3SJed Brown } 4079566063dSJacob Faibussowitsch PetscCall(PCGetOptionsPrefix(pc,&prefix)); 4089566063dSJacob Faibussowitsch PetscCall(PCGASMPrintSubdomains(pc)); 4098f3b4b4dSDmitry Karpeev 4106a4f0f73SDmitry Karpeev /* 4116a4f0f73SDmitry Karpeev Merge the ISs, create merged vectors and restrictions. 4126a4f0f73SDmitry Karpeev */ 4136a4f0f73SDmitry Karpeev /* Merge outer subdomain ISs and construct a restriction onto the disjoint union of local outer subdomains. */ 4146a4f0f73SDmitry Karpeev on = 0; 415f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 4169566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(osm->ois[i],&oni)); 4176a4f0f73SDmitry Karpeev on += oni; 418f746d493SDmitry Karpeev } 4199566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(on, &oidx)); 4206a4f0f73SDmitry Karpeev on = 0; 421930d09c1SFande Kong /* Merge local indices together */ 422f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 4239566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(osm->ois[i],&oni)); 4249566063dSJacob Faibussowitsch PetscCall(ISGetIndices(osm->ois[i],&oidxi)); 4259566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(oidx+on,oidxi,oni)); 4269566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(osm->ois[i],&oidxi)); 4276a4f0f73SDmitry Karpeev on += oni; 428f746d493SDmitry Karpeev } 4299566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(((PetscObject)(pc))->comm,on,oidx,PETSC_OWN_POINTER,&gois)); 430ea91fabdSFande Kong nTotalInnerIndices = 0; 431ea91fabdSFande Kong for (i=0; i<osm->n; i++) { 4329566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(osm->iis[i],&nInnerIndices)); 433ea91fabdSFande Kong nTotalInnerIndices += nInnerIndices; 434ea91fabdSFande Kong } 4359566063dSJacob Faibussowitsch PetscCall(VecCreateMPI(((PetscObject)(pc))->comm,nTotalInnerIndices,PETSC_DETERMINE,&x)); 4369566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x,&y)); 437ea91fabdSFande Kong 4389566063dSJacob Faibussowitsch PetscCall(VecCreateMPI(PetscObjectComm((PetscObject)pc),on,PETSC_DECIDE,&osm->gx)); 4399566063dSJacob Faibussowitsch PetscCall(VecDuplicate(osm->gx,&osm->gy)); 4409566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(osm->gx, &gostart, NULL)); 4419566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PetscObjectComm((PetscObject)pc),on,gostart,1, &goid)); 442930d09c1SFande Kong /* gois might indices not on local */ 4439566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(x,gois,osm->gx,goid, &(osm->gorestriction))); 4449566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(osm->n,&numbering)); 4459566063dSJacob Faibussowitsch PetscCall(PetscObjectsListGetGlobalNumbering(PetscObjectComm((PetscObject)pc),osm->n,(PetscObject*)osm->ois,NULL,numbering)); 4469566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 4479566063dSJacob Faibussowitsch PetscCall(ISDestroy(&gois)); 4482fa5cd67SKarl Rupp 4496a4f0f73SDmitry Karpeev /* Merge inner subdomain ISs and construct a restriction onto the disjoint union of local inner subdomains. */ 4502fa5cd67SKarl Rupp { 4512fa5cd67SKarl Rupp PetscInt ini; /* Number of indices the i-th a local inner subdomain. */ 4528966356dSPierre Jolivet PetscInt in; /* Number of indices in the disjoint union of local inner subdomains. */ 4536a4f0f73SDmitry Karpeev PetscInt *iidx; /* Global indices in the merged local inner subdomain. */ 4546a4f0f73SDmitry Karpeev PetscInt *ioidx; /* Global indices of the disjoint union of inner subdomains within the disjoint union of outer subdomains. */ 4556a4f0f73SDmitry Karpeev IS giis; /* IS for the disjoint union of inner subdomains. */ 4566a4f0f73SDmitry Karpeev IS giois; /* IS for the disjoint union of inner subdomains within the disjoint union of outer subdomains. */ 457f771a274SFande Kong PetscScalar *array; 458f771a274SFande Kong const PetscInt *indices; 459f771a274SFande Kong PetscInt k; 4606a4f0f73SDmitry Karpeev on = 0; 461f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 4629566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(osm->ois[i],&oni)); 4636a4f0f73SDmitry Karpeev on += oni; 464f746d493SDmitry Karpeev } 4659566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(on, &iidx)); 4669566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(on, &ioidx)); 4679566063dSJacob Faibussowitsch PetscCall(VecGetArray(y,&array)); 468e12b4729SFande Kong /* set communicator id to determine where overlap is */ 469f771a274SFande Kong in = 0; 470f771a274SFande Kong for (i=0; i<osm->n; i++) { 4719566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(osm->iis[i],&ini)); 472f771a274SFande Kong for (k = 0; k < ini; ++k) { 473f771a274SFande Kong array[in+k] = numbering[i]; 474f771a274SFande Kong } 475f771a274SFande Kong in += ini; 476f771a274SFande Kong } 4779566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(y,&array)); 4789566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->gorestriction,y,osm->gy,INSERT_VALUES,SCATTER_FORWARD)); 4799566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->gorestriction,y,osm->gy,INSERT_VALUES,SCATTER_FORWARD)); 4809566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(osm->gy,&gostart, NULL)); 4819566063dSJacob Faibussowitsch PetscCall(VecGetArray(osm->gy,&array)); 482f771a274SFande Kong on = 0; 483f771a274SFande Kong in = 0; 484f771a274SFande Kong for (i=0; i<osm->n; i++) { 4859566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(osm->ois[i],&oni)); 4869566063dSJacob Faibussowitsch PetscCall(ISGetIndices(osm->ois[i],&indices)); 487f771a274SFande Kong for (k=0; k<oni; k++) { 488e12b4729SFande Kong /* skip overlapping indices to get inner domain */ 48943081b6cSDmitry Karpeev if (PetscRealPart(array[on+k]) != numbering[i]) continue; 490f771a274SFande Kong iidx[in] = indices[k]; 491f771a274SFande Kong ioidx[in++] = gostart+on+k; 492f771a274SFande Kong } 4939566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(osm->ois[i], &indices)); 494f771a274SFande Kong on += oni; 495f771a274SFande Kong } 4969566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(osm->gy,&array)); 4979566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)pc),in,iidx,PETSC_OWN_POINTER,&giis)); 4989566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)pc),in,ioidx,PETSC_OWN_POINTER,&giois)); 4999566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(y,giis,osm->gy,giois,&osm->girestriction)); 5009566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 5019566063dSJacob Faibussowitsch PetscCall(ISDestroy(&giis)); 5029566063dSJacob Faibussowitsch PetscCall(ISDestroy(&giois)); 503b1a0a8a3SJed Brown } 5049566063dSJacob Faibussowitsch PetscCall(ISDestroy(&goid)); 5059566063dSJacob Faibussowitsch PetscCall(PetscFree(numbering)); 5062fa5cd67SKarl Rupp 507f746d493SDmitry Karpeev /* Create the subdomain work vectors. */ 5089566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(osm->n,&osm->x)); 5099566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(osm->n,&osm->y)); 5109566063dSJacob Faibussowitsch PetscCall(VecGetArray(osm->gx, &gxarray)); 5119566063dSJacob Faibussowitsch PetscCall(VecGetArray(osm->gy, &gyarray)); 5126a4f0f73SDmitry Karpeev for (i=0, on=0; i<osm->n; ++i, on += oni) { 5136a4f0f73SDmitry Karpeev PetscInt oNi; 5149566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(osm->ois[i],&oni)); 515930d09c1SFande Kong /* on a sub communicator */ 5169566063dSJacob Faibussowitsch PetscCall(ISGetSize(osm->ois[i],&oNi)); 5179566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(((PetscObject)(osm->ois[i]))->comm,1,oni,oNi,gxarray+on,&osm->x[i])); 5189566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(((PetscObject)(osm->ois[i]))->comm,1,oni,oNi,gyarray+on,&osm->y[i])); 519b1a0a8a3SJed Brown } 5209566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(osm->gx, &gxarray)); 5219566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(osm->gy, &gyarray)); 5228f3b4b4dSDmitry Karpeev /* Create the subdomain solvers */ 5239566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(osm->n,&osm->ksp)); 5248f3b4b4dSDmitry Karpeev for (i=0; i<osm->n; i++) { 5258f3b4b4dSDmitry Karpeev char subprefix[PETSC_MAX_PATH_LEN+1]; 5269566063dSJacob Faibussowitsch PetscCall(KSPCreate(((PetscObject)(osm->ois[i]))->comm,&ksp)); 5279566063dSJacob Faibussowitsch PetscCall(KSPSetErrorIfNotConverged(ksp,pc->erroriffailure)); 5289566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)ksp)); 5299566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,1)); 5309566063dSJacob Faibussowitsch PetscCall(KSPSetType(ksp,KSPPREONLY)); 5319566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp,&subpc)); /* Why do we need this here? */ 5328f3b4b4dSDmitry Karpeev if (subdomain_dm) { 5339566063dSJacob Faibussowitsch PetscCall(KSPSetDM(ksp,subdomain_dm[i])); 5349566063dSJacob Faibussowitsch PetscCall(DMDestroy(subdomain_dm+i)); 5358f3b4b4dSDmitry Karpeev } 5369566063dSJacob Faibussowitsch PetscCall(PCGetOptionsPrefix(pc,&prefix)); 5379566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(ksp,prefix)); 5388f3b4b4dSDmitry Karpeev if (subdomain_names && subdomain_names[i]) { 5399566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(subprefix,PETSC_MAX_PATH_LEN,"sub_%s_",subdomain_names[i])); 5409566063dSJacob Faibussowitsch PetscCall(KSPAppendOptionsPrefix(ksp,subprefix)); 5419566063dSJacob Faibussowitsch PetscCall(PetscFree(subdomain_names[i])); 5428f3b4b4dSDmitry Karpeev } 5439566063dSJacob Faibussowitsch PetscCall(KSPAppendOptionsPrefix(ksp,"sub_")); 544b1a0a8a3SJed Brown osm->ksp[i] = ksp; 545b1a0a8a3SJed Brown } 5469566063dSJacob Faibussowitsch PetscCall(PetscFree(subdomain_dm)); 5479566063dSJacob Faibussowitsch PetscCall(PetscFree(subdomain_names)); 548b1a0a8a3SJed Brown scall = MAT_INITIAL_MATRIX; 5498f3b4b4dSDmitry Karpeev } else { /* if (pc->setupcalled) */ 550b1a0a8a3SJed Brown /* 5518f3b4b4dSDmitry Karpeev Destroy the submatrices from the previous iteration 552b1a0a8a3SJed Brown */ 553b1a0a8a3SJed Brown if (pc->flag == DIFFERENT_NONZERO_PATTERN) { 5549566063dSJacob Faibussowitsch PetscCall(MatDestroyMatrices(osm->n,&osm->pmat)); 555b1a0a8a3SJed Brown scall = MAT_INITIAL_MATRIX; 556b1a0a8a3SJed Brown } 557ea91fabdSFande Kong if (osm->permutationIS) { 5589566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(pc->pmat,osm->permutationIS,osm->permutationIS,scall,&osm->permutationP)); 5599566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)osm->permutationP)); 560ea91fabdSFande Kong osm->pcmat = pc->pmat; 561ea91fabdSFande Kong pc->pmat = osm->permutationP; 562b1a0a8a3SJed Brown } 563ea91fabdSFande Kong } 564ea91fabdSFande Kong 565b1a0a8a3SJed Brown /* 5662da392ccSBarry Smith Extract the submatrices. 567b1a0a8a3SJed Brown */ 56881a5c4aaSDmitry Karpeev if (size > 1) { 5699566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatricesMPI(pc->pmat,osm->n,osm->ois,osm->ois,scall,&osm->pmat)); 5702fa5cd67SKarl Rupp } else { 5719566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(pc->pmat,osm->n,osm->ois,osm->ois,scall,&osm->pmat)); 57281a5c4aaSDmitry Karpeev } 573b1a0a8a3SJed Brown if (scall == MAT_INITIAL_MATRIX) { 5749566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptionsPrefix((PetscObject)pc->pmat,&pprefix)); 575f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 5769566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)osm->pmat[i])); 5779566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)osm->pmat[i],pprefix)); 578b1a0a8a3SJed Brown } 579b1a0a8a3SJed Brown } 580b1a0a8a3SJed Brown 581b1a0a8a3SJed Brown /* Return control to the user so that the submatrices can be modified (e.g., to apply 582b1a0a8a3SJed Brown different boundary conditions for the submatrices than for the global problem) */ 5839566063dSJacob Faibussowitsch PetscCall(PCModifySubMatrices(pc,osm->n,osm->ois,osm->ois,osm->pmat,pc->modifysubmatricesP)); 584b1a0a8a3SJed Brown 585b1a0a8a3SJed Brown /* 5866a4f0f73SDmitry Karpeev Loop over submatrices putting them into local ksps 587b1a0a8a3SJed Brown */ 588f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 5899566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(osm->ksp[i],osm->pmat[i],osm->pmat[i])); 5909566063dSJacob Faibussowitsch PetscCall(KSPGetOptionsPrefix(osm->ksp[i],&prefix)); 5919566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(osm->pmat[i],prefix)); 592b1a0a8a3SJed Brown if (!pc->setupcalled) { 5939566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(osm->ksp[i])); 594b1a0a8a3SJed Brown } 595b1a0a8a3SJed Brown } 596ea91fabdSFande Kong if (osm->pcmat) { 5979566063dSJacob Faibussowitsch PetscCall(MatDestroy(&pc->pmat)); 598ea91fabdSFande Kong pc->pmat = osm->pcmat; 5990a545947SLisandro Dalcin osm->pcmat = NULL; 600ea91fabdSFande Kong } 601b1a0a8a3SJed Brown PetscFunctionReturn(0); 602b1a0a8a3SJed Brown } 603b1a0a8a3SJed Brown 604f746d493SDmitry Karpeev static PetscErrorCode PCSetUpOnBlocks_GASM(PC pc) 605b1a0a8a3SJed Brown { 606f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 607b1a0a8a3SJed Brown PetscInt i; 608b1a0a8a3SJed Brown 609b1a0a8a3SJed Brown PetscFunctionBegin; 610f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 6119566063dSJacob Faibussowitsch PetscCall(KSPSetUp(osm->ksp[i])); 612b1a0a8a3SJed Brown } 613b1a0a8a3SJed Brown PetscFunctionReturn(0); 614b1a0a8a3SJed Brown } 615b1a0a8a3SJed Brown 616ea91fabdSFande Kong static PetscErrorCode PCApply_GASM(PC pc,Vec xin,Vec yout) 617b1a0a8a3SJed Brown { 618f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 619f746d493SDmitry Karpeev PetscInt i; 620ea91fabdSFande Kong Vec x,y; 621b1a0a8a3SJed Brown ScatterMode forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE; 622b1a0a8a3SJed Brown 623b1a0a8a3SJed Brown PetscFunctionBegin; 624ea91fabdSFande Kong if (osm->pctoouter) { 6259566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->pctoouter,xin,osm->pcx,INSERT_VALUES,SCATTER_REVERSE)); 6269566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->pctoouter,xin,osm->pcx,INSERT_VALUES,SCATTER_REVERSE)); 627ea91fabdSFande Kong x = osm->pcx; 628ea91fabdSFande Kong y = osm->pcy; 629ea91fabdSFande Kong } else { 630ea91fabdSFande Kong x = xin; 631ea91fabdSFande Kong y = yout; 632ea91fabdSFande Kong } 633b1a0a8a3SJed Brown /* 63448e38a8aSPierre Jolivet support for limiting the restriction or interpolation only to the inner 635b1a0a8a3SJed Brown subdomain values (leaving the other values 0). 636b1a0a8a3SJed Brown */ 637f746d493SDmitry Karpeev if (!(osm->type & PC_GASM_RESTRICT)) { 638b1a0a8a3SJed Brown /* have to zero the work RHS since scatter may leave some slots empty */ 6399566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(osm->gx)); 6409566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->girestriction,x,osm->gx,INSERT_VALUES,forward)); 6412fa5cd67SKarl Rupp } else { 6429566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward)); 643b1a0a8a3SJed Brown } 6449566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(osm->gy)); 6456a4f0f73SDmitry Karpeev if (!(osm->type & PC_GASM_RESTRICT)) { 6469566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->girestriction,x,osm->gx,INSERT_VALUES,forward)); 6472fa5cd67SKarl Rupp } else { 6489566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward)); 6496a4f0f73SDmitry Karpeev } 65086b96d36SDmitry Karpeev /* do the subdomain solves */ 65186b96d36SDmitry Karpeev for (i=0; i<osm->n; ++i) { 6529566063dSJacob Faibussowitsch PetscCall(KSPSolve(osm->ksp[i],osm->x[i],osm->y[i])); 6539566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(osm->ksp[i],pc,osm->y[i])); 654b1a0a8a3SJed Brown } 65548e38a8aSPierre Jolivet /* do we need to zero y? */ 6569566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(y)); 6576a4f0f73SDmitry Karpeev if (!(osm->type & PC_GASM_INTERPOLATE)) { 6589566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->girestriction,osm->gy,y,ADD_VALUES,reverse)); 6599566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->girestriction,osm->gy,y,ADD_VALUES,reverse)); 6602fa5cd67SKarl Rupp } else { 6619566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse)); 6629566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse)); 6636a4f0f73SDmitry Karpeev } 664ea91fabdSFande Kong if (osm->pctoouter) { 6659566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->pctoouter,y,yout,INSERT_VALUES,SCATTER_FORWARD)); 6669566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->pctoouter,y,yout,INSERT_VALUES,SCATTER_FORWARD)); 667ea91fabdSFande Kong } 668ea91fabdSFande Kong PetscFunctionReturn(0); 669b1a0a8a3SJed Brown } 670b1a0a8a3SJed Brown 67148e38a8aSPierre Jolivet static PetscErrorCode PCMatApply_GASM(PC pc,Mat Xin,Mat Yout) 67248e38a8aSPierre Jolivet { 67348e38a8aSPierre Jolivet PC_GASM *osm = (PC_GASM*)pc->data; 67448e38a8aSPierre Jolivet Mat X,Y,O=NULL,Z,W; 67548e38a8aSPierre Jolivet Vec x,y; 67648e38a8aSPierre Jolivet PetscInt i,m,M,N; 67748e38a8aSPierre Jolivet ScatterMode forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE; 67848e38a8aSPierre Jolivet 67948e38a8aSPierre Jolivet PetscFunctionBegin; 68008401ef6SPierre Jolivet PetscCheck(osm->n == 1,PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Not yet implemented"); 6819566063dSJacob Faibussowitsch PetscCall(MatGetSize(Xin,NULL,&N)); 68248e38a8aSPierre Jolivet if (osm->pctoouter) { 6839566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(osm->pcx,&m)); 6849566063dSJacob Faibussowitsch PetscCall(VecGetSize(osm->pcx,&M)); 6859566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PetscObjectComm((PetscObject)osm->ois[0]),m,PETSC_DECIDE,M,N,NULL,&O)); 68648e38a8aSPierre Jolivet for (i = 0; i < N; ++i) { 6879566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumnVecRead(Xin,i,&x)); 6889566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumnVecWrite(O,i,&y)); 6899566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->pctoouter,x,y,INSERT_VALUES,SCATTER_REVERSE)); 6909566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->pctoouter,x,y,INSERT_VALUES,SCATTER_REVERSE)); 6919566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumnVecWrite(O,i,&y)); 6929566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumnVecRead(Xin,i,&x)); 69348e38a8aSPierre Jolivet } 69448e38a8aSPierre Jolivet X = Y = O; 69548e38a8aSPierre Jolivet } else { 69648e38a8aSPierre Jolivet X = Xin; 69748e38a8aSPierre Jolivet Y = Yout; 69848e38a8aSPierre Jolivet } 69948e38a8aSPierre Jolivet /* 70048e38a8aSPierre Jolivet support for limiting the restriction or interpolation only to the inner 70148e38a8aSPierre Jolivet subdomain values (leaving the other values 0). 70248e38a8aSPierre Jolivet */ 7039566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(osm->x[0],&m)); 7049566063dSJacob Faibussowitsch PetscCall(VecGetSize(osm->x[0],&M)); 7059566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PetscObjectComm((PetscObject)osm->ois[0]),m,PETSC_DECIDE,M,N,NULL,&Z)); 70648e38a8aSPierre Jolivet for (i = 0; i < N; ++i) { 7079566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumnVecRead(X,i,&x)); 7089566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumnVecWrite(Z,i,&y)); 70948e38a8aSPierre Jolivet if (!(osm->type & PC_GASM_RESTRICT)) { 71048e38a8aSPierre Jolivet /* have to zero the work RHS since scatter may leave some slots empty */ 7119566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(y)); 7129566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->girestriction,x,y,INSERT_VALUES,forward)); 7139566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->girestriction,x,y,INSERT_VALUES,forward)); 71448e38a8aSPierre Jolivet } else { 7159566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->gorestriction,x,y,INSERT_VALUES,forward)); 7169566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->gorestriction,x,y,INSERT_VALUES,forward)); 71748e38a8aSPierre Jolivet } 7189566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumnVecWrite(Z,i,&y)); 7199566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumnVecRead(X,i,&x)); 72048e38a8aSPierre Jolivet } 7219566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PetscObjectComm((PetscObject)osm->ois[0]),m,PETSC_DECIDE,M,N,NULL,&W)); 7229566063dSJacob Faibussowitsch PetscCall(MatSetOption(Z,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE)); 7239566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(Z,MAT_FINAL_ASSEMBLY)); 7249566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(Z,MAT_FINAL_ASSEMBLY)); 72548e38a8aSPierre Jolivet /* do the subdomain solve */ 7269566063dSJacob Faibussowitsch PetscCall(KSPMatSolve(osm->ksp[0],Z,W)); 7279566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(osm->ksp[0],pc,NULL)); 7289566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Z)); 72948e38a8aSPierre Jolivet /* do we need to zero y? */ 7309566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(Y)); 73148e38a8aSPierre Jolivet for (i = 0; i < N; ++i) { 7329566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumnVecWrite(Y,i,&y)); 7339566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumnVecRead(W,i,&x)); 73448e38a8aSPierre Jolivet if (!(osm->type & PC_GASM_INTERPOLATE)) { 7359566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->girestriction,x,y,ADD_VALUES,reverse)); 7369566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->girestriction,x,y,ADD_VALUES,reverse)); 73748e38a8aSPierre Jolivet } else { 7389566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->gorestriction,x,y,ADD_VALUES,reverse)); 7399566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->gorestriction,x,y,ADD_VALUES,reverse)); 74048e38a8aSPierre Jolivet } 7419566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumnVecRead(W,i,&x)); 74248e38a8aSPierre Jolivet if (osm->pctoouter) { 7439566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumnVecWrite(Yout,i,&x)); 7449566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->pctoouter,y,x,INSERT_VALUES,SCATTER_FORWARD)); 7459566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->pctoouter,y,x,INSERT_VALUES,SCATTER_FORWARD)); 7469566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumnVecRead(Yout,i,&x)); 74748e38a8aSPierre Jolivet } 7489566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumnVecWrite(Y,i,&y)); 74948e38a8aSPierre Jolivet } 7509566063dSJacob Faibussowitsch PetscCall(MatDestroy(&W)); 7519566063dSJacob Faibussowitsch PetscCall(MatDestroy(&O)); 75248e38a8aSPierre Jolivet PetscFunctionReturn(0); 75348e38a8aSPierre Jolivet } 75448e38a8aSPierre Jolivet 755ea91fabdSFande Kong static PetscErrorCode PCApplyTranspose_GASM(PC pc,Vec xin,Vec yout) 756b1a0a8a3SJed Brown { 757f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 758f746d493SDmitry Karpeev PetscInt i; 759ea91fabdSFande Kong Vec x,y; 760b1a0a8a3SJed Brown ScatterMode forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE; 761b1a0a8a3SJed Brown 762b1a0a8a3SJed Brown PetscFunctionBegin; 763ea91fabdSFande Kong if (osm->pctoouter) { 7649566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->pctoouter,xin,osm->pcx,INSERT_VALUES,SCATTER_REVERSE)); 7659566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->pctoouter,xin,osm->pcx,INSERT_VALUES,SCATTER_REVERSE)); 766ea91fabdSFande Kong x = osm->pcx; 767ea91fabdSFande Kong y = osm->pcy; 768ea91fabdSFande Kong }else{ 769ea91fabdSFande Kong x = xin; 770ea91fabdSFande Kong y = yout; 771ea91fabdSFande Kong } 772b1a0a8a3SJed Brown /* 773b1a0a8a3SJed Brown Support for limiting the restriction or interpolation to only local 774b1a0a8a3SJed Brown subdomain values (leaving the other values 0). 775b1a0a8a3SJed Brown 776f746d493SDmitry Karpeev Note: these are reversed from the PCApply_GASM() because we are applying the 777b1a0a8a3SJed Brown transpose of the three terms 778b1a0a8a3SJed Brown */ 779f746d493SDmitry Karpeev if (!(osm->type & PC_GASM_INTERPOLATE)) { 780b1a0a8a3SJed Brown /* have to zero the work RHS since scatter may leave some slots empty */ 7819566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(osm->gx)); 7829566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->girestriction,x,osm->gx,INSERT_VALUES,forward)); 7832fa5cd67SKarl Rupp } else { 7849566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward)); 785b1a0a8a3SJed Brown } 7869566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(osm->gy)); 7876a4f0f73SDmitry Karpeev if (!(osm->type & PC_GASM_INTERPOLATE)) { 7889566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->girestriction,x,osm->gx,INSERT_VALUES,forward)); 7892fa5cd67SKarl Rupp } else { 7909566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward)); 7916a4f0f73SDmitry Karpeev } 792b1a0a8a3SJed Brown /* do the local solves */ 793ab3e923fSDmitry 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. */ 7949566063dSJacob Faibussowitsch PetscCall(KSPSolveTranspose(osm->ksp[i],osm->x[i],osm->y[i])); 7959566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(osm->ksp[i],pc,osm->y[i])); 796b1a0a8a3SJed Brown } 7979566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(y)); 7986a4f0f73SDmitry Karpeev if (!(osm->type & PC_GASM_RESTRICT)) { 7999566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->girestriction,osm->gy,y,ADD_VALUES,reverse)); 8009566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->girestriction,osm->gy,y,ADD_VALUES,reverse)); 8012fa5cd67SKarl Rupp } else { 8029566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse)); 8039566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse)); 8046a4f0f73SDmitry Karpeev } 805ea91fabdSFande Kong if (osm->pctoouter) { 8069566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->pctoouter,y,yout,INSERT_VALUES,SCATTER_FORWARD)); 8079566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->pctoouter,y,yout,INSERT_VALUES,SCATTER_FORWARD)); 808ea91fabdSFande Kong } 809b1a0a8a3SJed Brown PetscFunctionReturn(0); 810b1a0a8a3SJed Brown } 811b1a0a8a3SJed Brown 812a06653b4SBarry Smith static PetscErrorCode PCReset_GASM(PC pc) 813b1a0a8a3SJed Brown { 814f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 815b1a0a8a3SJed Brown PetscInt i; 816b1a0a8a3SJed Brown 817b1a0a8a3SJed Brown PetscFunctionBegin; 818b1a0a8a3SJed Brown if (osm->ksp) { 819f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 8209566063dSJacob Faibussowitsch PetscCall(KSPReset(osm->ksp[i])); 821b1a0a8a3SJed Brown } 822b1a0a8a3SJed Brown } 823b1a0a8a3SJed Brown if (osm->pmat) { 824f746d493SDmitry Karpeev if (osm->n > 0) { 825df750dc8SHong Zhang PetscMPIInt size; 8269566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size)); 827df750dc8SHong Zhang if (size > 1) { 8287dae84e0SHong Zhang /* osm->pmat is created by MatCreateSubMatricesMPI(), cannot use MatDestroySubMatrices() */ 8299566063dSJacob Faibussowitsch PetscCall(MatDestroyMatrices(osm->n,&osm->pmat)); 830df750dc8SHong Zhang } else { 8319566063dSJacob Faibussowitsch PetscCall(MatDestroySubMatrices(osm->n,&osm->pmat)); 832df750dc8SHong Zhang } 833b1a0a8a3SJed Brown } 834b1a0a8a3SJed Brown } 835d34fcf5fSBarry Smith if (osm->x) { 836f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 8379566063dSJacob Faibussowitsch PetscCall(VecDestroy(&osm->x[i])); 8389566063dSJacob Faibussowitsch PetscCall(VecDestroy(&osm->y[i])); 839b1a0a8a3SJed Brown } 840d34fcf5fSBarry Smith } 8419566063dSJacob Faibussowitsch PetscCall(VecDestroy(&osm->gx)); 8429566063dSJacob Faibussowitsch PetscCall(VecDestroy(&osm->gy)); 843ab3e923fSDmitry Karpeev 8449566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&osm->gorestriction)); 8459566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&osm->girestriction)); 8468f3b4b4dSDmitry Karpeev if (!osm->user_subdomains) { 8479566063dSJacob Faibussowitsch PetscCall(PCGASMDestroySubdomains(osm->n,&osm->ois,&osm->iis)); 8488f3b4b4dSDmitry Karpeev osm->N = PETSC_DETERMINE; 8498f3b4b4dSDmitry Karpeev osm->nmax = PETSC_DETERMINE; 8508f3b4b4dSDmitry Karpeev } 851ea91fabdSFande Kong if (osm->pctoouter) { 8529566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&(osm->pctoouter))); 853ea91fabdSFande Kong } 854ea91fabdSFande Kong if (osm->permutationIS) { 8559566063dSJacob Faibussowitsch PetscCall(ISDestroy(&(osm->permutationIS))); 856ea91fabdSFande Kong } 857ea91fabdSFande Kong if (osm->pcx) { 8589566063dSJacob Faibussowitsch PetscCall(VecDestroy(&(osm->pcx))); 859ea91fabdSFande Kong } 860ea91fabdSFande Kong if (osm->pcy) { 8619566063dSJacob Faibussowitsch PetscCall(VecDestroy(&(osm->pcy))); 862ea91fabdSFande Kong } 863ea91fabdSFande Kong if (osm->permutationP) { 8649566063dSJacob Faibussowitsch PetscCall(MatDestroy(&(osm->permutationP))); 865ea91fabdSFande Kong } 866ea91fabdSFande Kong if (osm->pcmat) { 8679566063dSJacob Faibussowitsch PetscCall(MatDestroy(&osm->pcmat)); 868ea91fabdSFande Kong } 869a06653b4SBarry Smith PetscFunctionReturn(0); 870a06653b4SBarry Smith } 871a06653b4SBarry Smith 872a06653b4SBarry Smith static PetscErrorCode PCDestroy_GASM(PC pc) 873a06653b4SBarry Smith { 874a06653b4SBarry Smith PC_GASM *osm = (PC_GASM*)pc->data; 875a06653b4SBarry Smith PetscInt i; 876a06653b4SBarry Smith 877a06653b4SBarry Smith PetscFunctionBegin; 8789566063dSJacob Faibussowitsch PetscCall(PCReset_GASM(pc)); 879135757f6SDmitry Karpeev /* PCReset will not destroy subdomains, if user_subdomains is true. */ 8809566063dSJacob Faibussowitsch PetscCall(PCGASMDestroySubdomains(osm->n,&osm->ois,&osm->iis)); 881a06653b4SBarry Smith if (osm->ksp) { 882a06653b4SBarry Smith for (i=0; i<osm->n; i++) { 8839566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&osm->ksp[i])); 884a06653b4SBarry Smith } 8859566063dSJacob Faibussowitsch PetscCall(PetscFree(osm->ksp)); 886a06653b4SBarry Smith } 8879566063dSJacob Faibussowitsch PetscCall(PetscFree(osm->x)); 8889566063dSJacob Faibussowitsch PetscCall(PetscFree(osm->y)); 8892e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetSubdomains_C",NULL)); 8902e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetOverlap_C",NULL)); 8912e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetType_C",NULL)); 8922e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetSortIndices_C",NULL)); 8932e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGASMGetSubKSP_C",NULL)); 8949566063dSJacob Faibussowitsch PetscCall(PetscFree(pc->data)); 895b1a0a8a3SJed Brown PetscFunctionReturn(0); 896b1a0a8a3SJed Brown } 897b1a0a8a3SJed Brown 8984416b707SBarry Smith static PetscErrorCode PCSetFromOptions_GASM(PetscOptionItems *PetscOptionsObject,PC pc) 899a6dfd86eSKarl Rupp { 900f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 901b1a0a8a3SJed Brown PetscInt blocks,ovl; 90257501b6eSBarry Smith PetscBool flg; 903f746d493SDmitry Karpeev PCGASMType gasmtype; 904b1a0a8a3SJed Brown 905b1a0a8a3SJed Brown PetscFunctionBegin; 906d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject,"Generalized additive Schwarz options"); 9079566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_gasm_use_dm_subdomains","If subdomains aren't set, use DMCreateDomainDecomposition() to define subdomains.","PCGASMSetUseDMSubdomains",osm->dm_subdomains,&osm->dm_subdomains,&flg)); 9089566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_gasm_total_subdomains","Total number of subdomains across communicator","PCGASMSetTotalSubdomains",osm->N,&blocks,&flg)); 909*1baa6e33SBarry Smith if (flg) PetscCall(PCGASMSetTotalSubdomains(pc,blocks)); 9109566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_gasm_overlap","Number of overlapping degrees of freedom","PCGASMSetOverlap",osm->overlap,&ovl,&flg)); 91165db9045SDmitry Karpeev if (flg) { 9129566063dSJacob Faibussowitsch PetscCall(PCGASMSetOverlap(pc,ovl)); 913d709ea83SDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 91465db9045SDmitry Karpeev } 915b1a0a8a3SJed Brown flg = PETSC_FALSE; 9169566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-pc_gasm_type","Type of restriction/extension","PCGASMSetType",PCGASMTypes,(PetscEnum)osm->type,(PetscEnum*)&gasmtype,&flg)); 9179566063dSJacob Faibussowitsch if (flg) PetscCall(PCGASMSetType(pc,gasmtype)); 9189566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_gasm_use_hierachical_partitioning","use hierarchical partitioning",NULL,osm->hierarchicalpartitioning,&osm->hierarchicalpartitioning,&flg)); 919d0609cedSBarry Smith PetscOptionsHeadEnd(); 920b1a0a8a3SJed Brown PetscFunctionReturn(0); 921b1a0a8a3SJed Brown } 922b1a0a8a3SJed Brown 923b1a0a8a3SJed Brown /*------------------------------------------------------------------------------------*/ 924b1a0a8a3SJed Brown 9258f3b4b4dSDmitry Karpeev /*@ 9268f3b4b4dSDmitry Karpeev PCGASMSetTotalSubdomains - sets the total number of subdomains to use across the 9278f3b4b4dSDmitry Karpeev communicator. 9288f3b4b4dSDmitry Karpeev Logically collective on pc 9298f3b4b4dSDmitry Karpeev 9308f3b4b4dSDmitry Karpeev Input Parameters: 9318f3b4b4dSDmitry Karpeev + pc - the preconditioner 9328f3b4b4dSDmitry Karpeev - N - total number of subdomains 9338f3b4b4dSDmitry Karpeev 9348f3b4b4dSDmitry Karpeev Level: beginner 9358f3b4b4dSDmitry Karpeev 936db781477SPatrick Sanan .seealso: `PCGASMSetSubdomains()`, `PCGASMSetOverlap()` 937db781477SPatrick Sanan `PCGASMCreateSubdomains2D()` 9388f3b4b4dSDmitry Karpeev @*/ 9398f3b4b4dSDmitry Karpeev PetscErrorCode PCGASMSetTotalSubdomains(PC pc,PetscInt N) 9408f3b4b4dSDmitry Karpeev { 9418f3b4b4dSDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 9428f3b4b4dSDmitry Karpeev PetscMPIInt size,rank; 9438f3b4b4dSDmitry Karpeev 9448f3b4b4dSDmitry Karpeev PetscFunctionBegin; 94563a3b9bcSJacob Faibussowitsch PetscCheck(N >= 1,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Total number of subdomains must be 1 or more, got N = %" PetscInt_FMT,N); 94628b400f6SJacob Faibussowitsch PetscCheck(!pc->setupcalled,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetTotalSubdomains() should be called before calling PCSetUp()."); 9478f3b4b4dSDmitry Karpeev 9489566063dSJacob Faibussowitsch PetscCall(PCGASMDestroySubdomains(osm->n,&osm->iis,&osm->ois)); 9498f3b4b4dSDmitry Karpeev osm->ois = osm->iis = NULL; 9508f3b4b4dSDmitry Karpeev 9519566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size)); 9529566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank)); 9538f3b4b4dSDmitry Karpeev osm->N = N; 9548f3b4b4dSDmitry Karpeev osm->n = PETSC_DETERMINE; 9558f3b4b4dSDmitry Karpeev osm->nmax = PETSC_DETERMINE; 9568f3b4b4dSDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 9578f3b4b4dSDmitry Karpeev PetscFunctionReturn(0); 9588f3b4b4dSDmitry Karpeev } 9598f3b4b4dSDmitry Karpeev 960b541e6a4SDmitry Karpeev static PetscErrorCode PCGASMSetSubdomains_GASM(PC pc,PetscInt n,IS iis[],IS ois[]) 961b1a0a8a3SJed Brown { 962f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 963b1a0a8a3SJed Brown PetscInt i; 964b1a0a8a3SJed Brown 965b1a0a8a3SJed Brown PetscFunctionBegin; 96663a3b9bcSJacob Faibussowitsch PetscCheck(n >= 1,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Each process must have 1 or more subdomains, got n = %" PetscInt_FMT,n); 96728b400f6SJacob Faibussowitsch PetscCheck(!pc->setupcalled,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetSubdomains() should be called before calling PCSetUp()."); 968b1a0a8a3SJed Brown 9699566063dSJacob Faibussowitsch PetscCall(PCGASMDestroySubdomains(osm->n,&osm->iis,&osm->ois)); 9708f3b4b4dSDmitry Karpeev osm->iis = osm->ois = NULL; 9718f3b4b4dSDmitry Karpeev osm->n = n; 9728f3b4b4dSDmitry Karpeev osm->N = PETSC_DETERMINE; 9738f3b4b4dSDmitry Karpeev osm->nmax = PETSC_DETERMINE; 974a35b7b57SDmitry Karpeev if (ois) { 9759566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n,&osm->ois)); 9768f3b4b4dSDmitry Karpeev for (i=0; i<n; i++) { 9779566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)ois[i])); 9788f3b4b4dSDmitry Karpeev osm->ois[i] = ois[i]; 9798f3b4b4dSDmitry Karpeev } 9808f3b4b4dSDmitry Karpeev /* 9818f3b4b4dSDmitry Karpeev Since the user set the outer subdomains, even if nontrivial overlap was requested via PCGASMSetOverlap(), 9828f3b4b4dSDmitry Karpeev it will be ignored. To avoid confusion later on (e.g., when viewing the PC), the overlap size is set to -1. 9838f3b4b4dSDmitry Karpeev */ 984b1a0a8a3SJed Brown osm->overlap = -1; 985670417b2SFande Kong /* inner subdomains must be provided */ 98628b400f6SJacob Faibussowitsch PetscCheck(iis,PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"inner indices have to be provided "); 987670417b2SFande Kong }/* end if */ 988a35b7b57SDmitry Karpeev if (iis) { 9899566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n,&osm->iis)); 9908f3b4b4dSDmitry Karpeev for (i=0; i<n; i++) { 9919566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)iis[i])); 9928f3b4b4dSDmitry Karpeev osm->iis[i] = iis[i]; 9938f3b4b4dSDmitry Karpeev } 994a35b7b57SDmitry Karpeev if (!ois) { 995390e1bf2SBarry Smith osm->ois = NULL; 996670417b2SFande Kong /* if user does not provide outer indices, we will create the corresponding outer indices using osm->overlap =1 in PCSetUp_GASM */ 997670417b2SFande Kong } 998670417b2SFande Kong } 99976bd3646SJed Brown if (PetscDefined(USE_DEBUG)) { 1000670417b2SFande Kong PetscInt j,rstart,rend,*covered,lsize; 1001670417b2SFande Kong const PetscInt *indices; 1002670417b2SFande Kong /* check if the inner indices cover and only cover the local portion of the preconditioning matrix */ 10039566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(pc->pmat,&rstart,&rend)); 10049566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(rend-rstart,&covered)); 1005670417b2SFande Kong /* check if the current processor owns indices from others */ 1006a35b7b57SDmitry Karpeev for (i=0; i<n; i++) { 10079566063dSJacob Faibussowitsch PetscCall(ISGetIndices(osm->iis[i],&indices)); 10089566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(osm->iis[i],&lsize)); 1009670417b2SFande Kong for (j=0; j<lsize; j++) { 10102472a847SBarry Smith PetscCheck(indices[j] >= rstart && indices[j] < rend,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"inner subdomains can not own an index %" PetscInt_FMT " from other processors", indices[j]); 10112472a847SBarry Smith PetscCheck(covered[indices[j]-rstart] != 1,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"inner subdomains can not have an overlapping index %" PetscInt_FMT " ",indices[j]); 101263a3b9bcSJacob Faibussowitsch covered[indices[j]-rstart] = 1; 1013a35b7b57SDmitry Karpeev } 10149566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(osm->iis[i],&indices)); 10158f3b4b4dSDmitry Karpeev } 1016670417b2SFande Kong /* check if we miss any indices */ 1017670417b2SFande Kong for (i=rstart; i<rend; i++) { 101863a3b9bcSJacob Faibussowitsch PetscCheck(covered[i-rstart],PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"local entity %" PetscInt_FMT " was not covered by inner subdomains",i); 1019a35b7b57SDmitry Karpeev } 10209566063dSJacob Faibussowitsch PetscCall(PetscFree(covered)); 1021a35b7b57SDmitry Karpeev } 10228f3b4b4dSDmitry Karpeev if (iis) osm->user_subdomains = PETSC_TRUE; 1023b1a0a8a3SJed Brown PetscFunctionReturn(0); 1024b1a0a8a3SJed Brown } 1025b1a0a8a3SJed Brown 1026f7a08781SBarry Smith static PetscErrorCode PCGASMSetOverlap_GASM(PC pc,PetscInt ovl) 1027b1a0a8a3SJed Brown { 1028f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 1029b1a0a8a3SJed Brown 1030b1a0a8a3SJed Brown PetscFunctionBegin; 103108401ef6SPierre Jolivet PetscCheck(ovl >= 0,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap value requested"); 10322472a847SBarry Smith PetscCheck(!pc->setupcalled || ovl == osm->overlap,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetOverlap() should be called before PCSetUp()."); 10332fa5cd67SKarl Rupp if (!pc->setupcalled) osm->overlap = ovl; 1034b1a0a8a3SJed Brown PetscFunctionReturn(0); 1035b1a0a8a3SJed Brown } 1036b1a0a8a3SJed Brown 1037f7a08781SBarry Smith static PetscErrorCode PCGASMSetType_GASM(PC pc,PCGASMType type) 1038b1a0a8a3SJed Brown { 1039f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 1040b1a0a8a3SJed Brown 1041b1a0a8a3SJed Brown PetscFunctionBegin; 1042b1a0a8a3SJed Brown osm->type = type; 1043b1a0a8a3SJed Brown osm->type_set = PETSC_TRUE; 1044b1a0a8a3SJed Brown PetscFunctionReturn(0); 1045b1a0a8a3SJed Brown } 1046b1a0a8a3SJed Brown 1047f7a08781SBarry Smith static PetscErrorCode PCGASMSetSortIndices_GASM(PC pc,PetscBool doSort) 1048b1a0a8a3SJed Brown { 1049f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 1050b1a0a8a3SJed Brown 1051b1a0a8a3SJed Brown PetscFunctionBegin; 1052b1a0a8a3SJed Brown osm->sort_indices = doSort; 1053b1a0a8a3SJed Brown PetscFunctionReturn(0); 1054b1a0a8a3SJed Brown } 1055b1a0a8a3SJed Brown 105644fe18e5SDmitry Karpeev /* 10578f3b4b4dSDmitry Karpeev FIXME: This routine might need to be modified now that multiple ranks per subdomain are allowed. 105844fe18e5SDmitry Karpeev In particular, it would upset the global subdomain number calculation. 105944fe18e5SDmitry Karpeev */ 1060f7a08781SBarry Smith static PetscErrorCode PCGASMGetSubKSP_GASM(PC pc,PetscInt *n,PetscInt *first,KSP **ksp) 1061b1a0a8a3SJed Brown { 1062f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 1063b1a0a8a3SJed Brown 1064b1a0a8a3SJed Brown PetscFunctionBegin; 106508401ef6SPierre Jolivet PetscCheck(osm->n >= 1,PetscObjectComm((PetscObject)pc),PETSC_ERR_ORDER,"Need to call PCSetUp() on PC (or KSPSetUp() on the outer KSP object) before calling here"); 1066b1a0a8a3SJed Brown 10672fa5cd67SKarl Rupp if (n) *n = osm->n; 1068ab3e923fSDmitry Karpeev if (first) { 10699566063dSJacob Faibussowitsch PetscCallMPI(MPI_Scan(&osm->n,first,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)pc))); 1070ab3e923fSDmitry Karpeev *first -= osm->n; 1071b1a0a8a3SJed Brown } 1072b1a0a8a3SJed Brown if (ksp) { 1073b1a0a8a3SJed Brown /* Assume that local solves are now different; not necessarily 107406b43e7eSDmitry Karpeev true, though! This flag is used only for PCView_GASM() */ 1075b1a0a8a3SJed Brown *ksp = osm->ksp; 10766a4f0f73SDmitry Karpeev osm->same_subdomain_solvers = PETSC_FALSE; 1077b1a0a8a3SJed Brown } 1078b1a0a8a3SJed Brown PetscFunctionReturn(0); 1079ab3e923fSDmitry Karpeev } /* PCGASMGetSubKSP_GASM() */ 1080b1a0a8a3SJed Brown 1081b1a0a8a3SJed Brown /*@C 108206b43e7eSDmitry Karpeev PCGASMSetSubdomains - Sets the subdomains for this processor 108306b43e7eSDmitry Karpeev for the additive Schwarz preconditioner. 1084b1a0a8a3SJed Brown 10858f3b4b4dSDmitry Karpeev Collective on pc 1086b1a0a8a3SJed Brown 1087b1a0a8a3SJed Brown Input Parameters: 10888f3b4b4dSDmitry Karpeev + pc - the preconditioner object 108906b43e7eSDmitry Karpeev . n - the number of subdomains for this processor 10908f3b4b4dSDmitry Karpeev . iis - the index sets that define the inner subdomains (or NULL for PETSc to determine subdomains) 10918f3b4b4dSDmitry Karpeev - ois - the index sets that define the outer subdomains (or NULL to use the same as iis, or to construct by expanding iis by the requested overlap) 1092b1a0a8a3SJed Brown 1093b1a0a8a3SJed Brown Notes: 109406b43e7eSDmitry Karpeev The IS indices use the parallel, global numbering of the vector entries. 10956a4f0f73SDmitry Karpeev Inner subdomains are those where the correction is applied. 10966a4f0f73SDmitry Karpeev Outer subdomains are those where the residual necessary to obtain the 10976a4f0f73SDmitry Karpeev corrections is obtained (see PCGASMType for the use of inner/outer subdomains). 10986a4f0f73SDmitry Karpeev Both inner and outer subdomains can extend over several processors. 10996a4f0f73SDmitry Karpeev This processor's portion of a subdomain is known as a local subdomain. 11006a4f0f73SDmitry Karpeev 1101670417b2SFande Kong Inner subdomains can not overlap with each other, do not have any entities from remote processors, 1102670417b2SFande Kong and have to cover the entire local subdomain owned by the current processor. The index sets on each 1103670417b2SFande Kong process should be ordered such that the ith local subdomain is connected to the ith remote subdomain 1104670417b2SFande Kong on another MPI process. 1105670417b2SFande Kong 11066a4f0f73SDmitry Karpeev By default the GASM preconditioner uses 1 (local) subdomain per processor. 11076a4f0f73SDmitry Karpeev 1108b1a0a8a3SJed Brown Level: advanced 1109b1a0a8a3SJed Brown 1110db781477SPatrick Sanan .seealso: `PCGASMSetOverlap()`, `PCGASMGetSubKSP()`, 1111db781477SPatrick Sanan `PCGASMCreateSubdomains2D()`, `PCGASMGetSubdomains()` 1112b1a0a8a3SJed Brown @*/ 11136a4f0f73SDmitry Karpeev PetscErrorCode PCGASMSetSubdomains(PC pc,PetscInt n,IS iis[],IS ois[]) 1114b1a0a8a3SJed Brown { 11158f3b4b4dSDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 1116b1a0a8a3SJed Brown 1117b1a0a8a3SJed Brown PetscFunctionBegin; 1118b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1119cac4c232SBarry Smith PetscTryMethod(pc,"PCGASMSetSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,n,iis,ois)); 11208f3b4b4dSDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 1121b1a0a8a3SJed Brown PetscFunctionReturn(0); 1122b1a0a8a3SJed Brown } 1123b1a0a8a3SJed Brown 1124b1a0a8a3SJed Brown /*@ 1125f746d493SDmitry Karpeev PCGASMSetOverlap - Sets the overlap between a pair of subdomains for the 1126b1a0a8a3SJed Brown additive Schwarz preconditioner. Either all or no processors in the 11278f3b4b4dSDmitry Karpeev pc communicator must call this routine. 1128b1a0a8a3SJed Brown 11298f3b4b4dSDmitry Karpeev Logically Collective on pc 1130b1a0a8a3SJed Brown 1131b1a0a8a3SJed Brown Input Parameters: 1132b1a0a8a3SJed Brown + pc - the preconditioner context 11338f3b4b4dSDmitry Karpeev - ovl - the amount of overlap between subdomains (ovl >= 0, default value = 0) 1134b1a0a8a3SJed Brown 1135b1a0a8a3SJed Brown Options Database Key: 113606b43e7eSDmitry Karpeev . -pc_gasm_overlap <overlap> - Sets overlap 1137b1a0a8a3SJed Brown 1138b1a0a8a3SJed Brown Notes: 113906b43e7eSDmitry Karpeev By default the GASM preconditioner uses 1 subdomain per processor. To use 11408f3b4b4dSDmitry Karpeev multiple subdomain per perocessor or "straddling" subdomains that intersect 11418f3b4b4dSDmitry Karpeev multiple processors use PCGASMSetSubdomains() (or option -pc_gasm_total_subdomains <n>). 1142b1a0a8a3SJed Brown 11438f3b4b4dSDmitry Karpeev The overlap defaults to 0, so if one desires that no additional 1144b1a0a8a3SJed Brown overlap be computed beyond what may have been set with a call to 11458f3b4b4dSDmitry Karpeev PCGASMSetSubdomains(), then ovl must be set to be 0. In particular, if one does 11468f3b4b4dSDmitry Karpeev not explicitly set the subdomains in application code, then all overlap would be computed 1147f746d493SDmitry Karpeev internally by PETSc, and using an overlap of 0 would result in an GASM 1148b1a0a8a3SJed Brown variant that is equivalent to the block Jacobi preconditioner. 1149b1a0a8a3SJed Brown 1150b1a0a8a3SJed Brown Note that one can define initial index sets with any overlap via 115106b43e7eSDmitry Karpeev PCGASMSetSubdomains(); the routine PCGASMSetOverlap() merely allows 115206b43e7eSDmitry Karpeev PETSc to extend that overlap further, if desired. 1153b1a0a8a3SJed Brown 1154b1a0a8a3SJed Brown Level: intermediate 1155b1a0a8a3SJed Brown 1156db781477SPatrick Sanan .seealso: `PCGASMSetSubdomains()`, `PCGASMGetSubKSP()`, 1157db781477SPatrick Sanan `PCGASMCreateSubdomains2D()`, `PCGASMGetSubdomains()` 1158b1a0a8a3SJed Brown @*/ 11597087cfbeSBarry Smith PetscErrorCode PCGASMSetOverlap(PC pc,PetscInt ovl) 1160b1a0a8a3SJed Brown { 11618f3b4b4dSDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 1162b1a0a8a3SJed Brown 1163b1a0a8a3SJed Brown PetscFunctionBegin; 1164b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1165b1a0a8a3SJed Brown PetscValidLogicalCollectiveInt(pc,ovl,2); 1166cac4c232SBarry Smith PetscTryMethod(pc,"PCGASMSetOverlap_C",(PC,PetscInt),(pc,ovl)); 11678f3b4b4dSDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 1168b1a0a8a3SJed Brown PetscFunctionReturn(0); 1169b1a0a8a3SJed Brown } 1170b1a0a8a3SJed Brown 1171b1a0a8a3SJed Brown /*@ 1172f746d493SDmitry Karpeev PCGASMSetType - Sets the type of restriction and interpolation used 1173b1a0a8a3SJed Brown for local problems in the additive Schwarz method. 1174b1a0a8a3SJed Brown 1175b1a0a8a3SJed Brown Logically Collective on PC 1176b1a0a8a3SJed Brown 1177b1a0a8a3SJed Brown Input Parameters: 1178b1a0a8a3SJed Brown + pc - the preconditioner context 1179f746d493SDmitry Karpeev - type - variant of GASM, one of 1180b1a0a8a3SJed Brown .vb 1181f746d493SDmitry Karpeev PC_GASM_BASIC - full interpolation and restriction 1182f746d493SDmitry Karpeev PC_GASM_RESTRICT - full restriction, local processor interpolation 1183f746d493SDmitry Karpeev PC_GASM_INTERPOLATE - full interpolation, local processor restriction 1184f746d493SDmitry Karpeev PC_GASM_NONE - local processor restriction and interpolation 1185b1a0a8a3SJed Brown .ve 1186b1a0a8a3SJed Brown 1187b1a0a8a3SJed Brown Options Database Key: 1188f746d493SDmitry Karpeev . -pc_gasm_type [basic,restrict,interpolate,none] - Sets GASM type 1189b1a0a8a3SJed Brown 1190b1a0a8a3SJed Brown Level: intermediate 1191b1a0a8a3SJed Brown 1192db781477SPatrick Sanan .seealso: `PCGASMSetSubdomains()`, `PCGASMGetSubKSP()`, 1193db781477SPatrick Sanan `PCGASMCreateSubdomains2D()` 1194b1a0a8a3SJed Brown @*/ 11957087cfbeSBarry Smith PetscErrorCode PCGASMSetType(PC pc,PCGASMType type) 1196b1a0a8a3SJed Brown { 1197b1a0a8a3SJed Brown PetscFunctionBegin; 1198b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1199b1a0a8a3SJed Brown PetscValidLogicalCollectiveEnum(pc,type,2); 1200cac4c232SBarry Smith PetscTryMethod(pc,"PCGASMSetType_C",(PC,PCGASMType),(pc,type)); 1201b1a0a8a3SJed Brown PetscFunctionReturn(0); 1202b1a0a8a3SJed Brown } 1203b1a0a8a3SJed Brown 1204b1a0a8a3SJed Brown /*@ 1205f746d493SDmitry Karpeev PCGASMSetSortIndices - Determines whether subdomain indices are sorted. 1206b1a0a8a3SJed Brown 1207b1a0a8a3SJed Brown Logically Collective on PC 1208b1a0a8a3SJed Brown 1209b1a0a8a3SJed Brown Input Parameters: 1210b1a0a8a3SJed Brown + pc - the preconditioner context 1211b1a0a8a3SJed Brown - doSort - sort the subdomain indices 1212b1a0a8a3SJed Brown 1213b1a0a8a3SJed Brown Level: intermediate 1214b1a0a8a3SJed Brown 1215db781477SPatrick Sanan .seealso: `PCGASMSetSubdomains()`, `PCGASMGetSubKSP()`, 1216db781477SPatrick Sanan `PCGASMCreateSubdomains2D()` 1217b1a0a8a3SJed Brown @*/ 12187087cfbeSBarry Smith PetscErrorCode PCGASMSetSortIndices(PC pc,PetscBool doSort) 1219b1a0a8a3SJed Brown { 1220b1a0a8a3SJed Brown PetscFunctionBegin; 1221b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1222b1a0a8a3SJed Brown PetscValidLogicalCollectiveBool(pc,doSort,2); 1223cac4c232SBarry Smith PetscTryMethod(pc,"PCGASMSetSortIndices_C",(PC,PetscBool),(pc,doSort)); 1224b1a0a8a3SJed Brown PetscFunctionReturn(0); 1225b1a0a8a3SJed Brown } 1226b1a0a8a3SJed Brown 1227b1a0a8a3SJed Brown /*@C 1228f746d493SDmitry Karpeev PCGASMGetSubKSP - Gets the local KSP contexts for all blocks on 1229b1a0a8a3SJed Brown this processor. 1230b1a0a8a3SJed Brown 1231b1a0a8a3SJed Brown Collective on PC iff first_local is requested 1232b1a0a8a3SJed Brown 1233b1a0a8a3SJed Brown Input Parameter: 1234b1a0a8a3SJed Brown . pc - the preconditioner context 1235b1a0a8a3SJed Brown 1236b1a0a8a3SJed Brown Output Parameters: 12370298fd71SBarry Smith + n_local - the number of blocks on this processor or NULL 12380298fd71SBarry Smith . first_local - the global number of the first block on this processor or NULL, 12390298fd71SBarry Smith all processors must request or all must pass NULL 1240b1a0a8a3SJed Brown - ksp - the array of KSP contexts 1241b1a0a8a3SJed Brown 1242b1a0a8a3SJed Brown Note: 1243f746d493SDmitry Karpeev After PCGASMGetSubKSP() the array of KSPes is not to be freed 1244b1a0a8a3SJed Brown 1245b1a0a8a3SJed Brown Currently for some matrix implementations only 1 block per processor 1246b1a0a8a3SJed Brown is supported. 1247b1a0a8a3SJed Brown 1248f746d493SDmitry Karpeev You must call KSPSetUp() before calling PCGASMGetSubKSP(). 1249b1a0a8a3SJed Brown 1250b1a0a8a3SJed Brown Level: advanced 1251b1a0a8a3SJed Brown 1252db781477SPatrick Sanan .seealso: `PCGASMSetSubdomains()`, `PCGASMSetOverlap()`, 1253db781477SPatrick Sanan `PCGASMCreateSubdomains2D()`, 1254b1a0a8a3SJed Brown @*/ 12557087cfbeSBarry Smith PetscErrorCode PCGASMGetSubKSP(PC pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp[]) 1256b1a0a8a3SJed Brown { 1257b1a0a8a3SJed Brown PetscFunctionBegin; 1258b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1259cac4c232SBarry Smith PetscUseMethod(pc,"PCGASMGetSubKSP_C",(PC,PetscInt*,PetscInt*,KSP **),(pc,n_local,first_local,ksp)); 1260b1a0a8a3SJed Brown PetscFunctionReturn(0); 1261b1a0a8a3SJed Brown } 1262b1a0a8a3SJed Brown 1263b1a0a8a3SJed Brown /* -------------------------------------------------------------------------------------*/ 1264b1a0a8a3SJed Brown /*MC 1265f746d493SDmitry Karpeev PCGASM - Use the (restricted) additive Schwarz method, each block is (approximately) solved with 1266b1a0a8a3SJed Brown its own KSP object. 1267b1a0a8a3SJed Brown 1268b1a0a8a3SJed Brown Options Database Keys: 12698f3b4b4dSDmitry Karpeev + -pc_gasm_total_subdomains <n> - Sets total number of local subdomains to be distributed among processors 127006b43e7eSDmitry Karpeev . -pc_gasm_view_subdomains - activates the printing of subdomain indices in PCView(), -ksp_view or -snes_view 127106b43e7eSDmitry Karpeev . -pc_gasm_print_subdomains - activates the printing of subdomain indices in PCSetUp() 127206b43e7eSDmitry Karpeev . -pc_gasm_overlap <ovl> - Sets overlap by which to (automatically) extend local subdomains 1273f746d493SDmitry Karpeev - -pc_gasm_type [basic,restrict,interpolate,none] - Sets GASM type 1274b1a0a8a3SJed Brown 1275b1a0a8a3SJed Brown IMPORTANT: If you run with, for example, 3 blocks on 1 processor or 3 blocks on 3 processors you 1276f746d493SDmitry Karpeev will get a different convergence rate due to the default option of -pc_gasm_type restrict. Use 1277f746d493SDmitry Karpeev -pc_gasm_type basic to use the standard GASM. 1278b1a0a8a3SJed Brown 127995452b02SPatrick Sanan Notes: 128095452b02SPatrick Sanan Blocks can be shared by multiple processes. 1281b1a0a8a3SJed Brown 1282b1a0a8a3SJed Brown To set options on the solvers for each block append -sub_ to all the KSP, and PC 1283b1a0a8a3SJed Brown options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly 1284b1a0a8a3SJed Brown 1285f746d493SDmitry Karpeev To set the options on the solvers separate for each block call PCGASMGetSubKSP() 1286b1a0a8a3SJed Brown and set the options directly on the resulting KSP object (you can access its PC 1287b1a0a8a3SJed Brown with KSPGetPC()) 1288b1a0a8a3SJed Brown 1289b1a0a8a3SJed Brown Level: beginner 1290b1a0a8a3SJed Brown 1291b1a0a8a3SJed Brown References: 1292606c0280SSatish Balay + * - M Dryja, OB Widlund, An additive variant of the Schwarz alternating method for the case of many subregions 129396a0c994SBarry Smith Courant Institute, New York University Technical report 1294606c0280SSatish Balay - * - Barry Smith, Petter Bjorstad, and William Gropp, Domain Decompositions: Parallel Multilevel Methods for Elliptic Partial Differential Equations, 129596a0c994SBarry Smith Cambridge University Press. 1296b1a0a8a3SJed Brown 1297db781477SPatrick Sanan .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, 1298db781477SPatrick Sanan `PCBJACOBI`, `PCGASMGetSubKSP()`, `PCGASMSetSubdomains()`, 1299db781477SPatrick Sanan `PCSetModifySubMatrices()`, `PCGASMSetOverlap()`, `PCGASMSetType()` 1300b1a0a8a3SJed Brown 1301b1a0a8a3SJed Brown M*/ 1302b1a0a8a3SJed Brown 13038cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_GASM(PC pc) 1304b1a0a8a3SJed Brown { 1305f746d493SDmitry Karpeev PC_GASM *osm; 1306b1a0a8a3SJed Brown 1307b1a0a8a3SJed Brown PetscFunctionBegin; 13089566063dSJacob Faibussowitsch PetscCall(PetscNewLog(pc,&osm)); 13092fa5cd67SKarl Rupp 13108f3b4b4dSDmitry Karpeev osm->N = PETSC_DETERMINE; 131106b43e7eSDmitry Karpeev osm->n = PETSC_DECIDE; 13128f3b4b4dSDmitry Karpeev osm->nmax = PETSC_DETERMINE; 13138f3b4b4dSDmitry Karpeev osm->overlap = 0; 13140a545947SLisandro Dalcin osm->ksp = NULL; 13150a545947SLisandro Dalcin osm->gorestriction = NULL; 13160a545947SLisandro Dalcin osm->girestriction = NULL; 13170a545947SLisandro Dalcin osm->pctoouter = NULL; 13180a545947SLisandro Dalcin osm->gx = NULL; 13190a545947SLisandro Dalcin osm->gy = NULL; 13200a545947SLisandro Dalcin osm->x = NULL; 13210a545947SLisandro Dalcin osm->y = NULL; 13220a545947SLisandro Dalcin osm->pcx = NULL; 13230a545947SLisandro Dalcin osm->pcy = NULL; 13240a545947SLisandro Dalcin osm->permutationIS = NULL; 13250a545947SLisandro Dalcin osm->permutationP = NULL; 13260a545947SLisandro Dalcin osm->pcmat = NULL; 13270a545947SLisandro Dalcin osm->ois = NULL; 13280a545947SLisandro Dalcin osm->iis = NULL; 13290a545947SLisandro Dalcin osm->pmat = NULL; 1330f746d493SDmitry Karpeev osm->type = PC_GASM_RESTRICT; 13316a4f0f73SDmitry Karpeev osm->same_subdomain_solvers = PETSC_TRUE; 1332b1a0a8a3SJed Brown osm->sort_indices = PETSC_TRUE; 1333d709ea83SDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 1334ea91fabdSFande Kong osm->hierarchicalpartitioning = PETSC_FALSE; 1335b1a0a8a3SJed Brown 1336b1a0a8a3SJed Brown pc->data = (void*)osm; 1337f746d493SDmitry Karpeev pc->ops->apply = PCApply_GASM; 133848e38a8aSPierre Jolivet pc->ops->matapply = PCMatApply_GASM; 1339f746d493SDmitry Karpeev pc->ops->applytranspose = PCApplyTranspose_GASM; 1340f746d493SDmitry Karpeev pc->ops->setup = PCSetUp_GASM; 1341a06653b4SBarry Smith pc->ops->reset = PCReset_GASM; 1342f746d493SDmitry Karpeev pc->ops->destroy = PCDestroy_GASM; 1343f746d493SDmitry Karpeev pc->ops->setfromoptions = PCSetFromOptions_GASM; 1344f746d493SDmitry Karpeev pc->ops->setuponblocks = PCSetUpOnBlocks_GASM; 1345f746d493SDmitry Karpeev pc->ops->view = PCView_GASM; 13460a545947SLisandro Dalcin pc->ops->applyrichardson = NULL; 1347b1a0a8a3SJed Brown 13489566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetSubdomains_C",PCGASMSetSubdomains_GASM)); 13499566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetOverlap_C",PCGASMSetOverlap_GASM)); 13509566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetType_C",PCGASMSetType_GASM)); 13519566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetSortIndices_C",PCGASMSetSortIndices_GASM)); 13529566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGASMGetSubKSP_C",PCGASMGetSubKSP_GASM)); 1353b1a0a8a3SJed Brown PetscFunctionReturn(0); 1354b1a0a8a3SJed Brown } 1355b1a0a8a3SJed Brown 13568f3b4b4dSDmitry Karpeev PetscErrorCode PCGASMCreateLocalSubdomains(Mat A, PetscInt nloc, IS *iis[]) 1357b1a0a8a3SJed Brown { 1358b1a0a8a3SJed Brown MatPartitioning mpart; 1359b1a0a8a3SJed Brown const char *prefix; 1360b1a0a8a3SJed Brown PetscInt i,j,rstart,rend,bs; 1361976e8c5aSLisandro Dalcin PetscBool hasop, isbaij = PETSC_FALSE,foundpart = PETSC_FALSE; 13620298fd71SBarry Smith Mat Ad = NULL, adj; 1363b1a0a8a3SJed Brown IS ispart,isnumb,*is; 1364b1a0a8a3SJed Brown 1365b1a0a8a3SJed Brown PetscFunctionBegin; 136663a3b9bcSJacob Faibussowitsch PetscCheck(nloc >= 1,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of local subdomains must > 0, got nloc = %" PetscInt_FMT,nloc); 1367b1a0a8a3SJed Brown 1368b1a0a8a3SJed Brown /* Get prefix, row distribution, and block size */ 13699566063dSJacob Faibussowitsch PetscCall(MatGetOptionsPrefix(A,&prefix)); 13709566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A,&rstart,&rend)); 13719566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(A,&bs)); 13722472a847SBarry Smith PetscCheck(rstart/bs*bs == rstart && rend/bs*bs == rend,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"bad row distribution [%" PetscInt_FMT ",%" PetscInt_FMT ") for matrix block size %" PetscInt_FMT,rstart,rend,bs); 1373b1a0a8a3SJed Brown 1374b1a0a8a3SJed Brown /* Get diagonal block from matrix if possible */ 13759566063dSJacob Faibussowitsch PetscCall(MatHasOperation(A,MATOP_GET_DIAGONAL_BLOCK,&hasop)); 1376976e8c5aSLisandro Dalcin if (hasop) { 13779566063dSJacob Faibussowitsch PetscCall(MatGetDiagonalBlock(A,&Ad)); 1378b1a0a8a3SJed Brown } 1379b1a0a8a3SJed Brown if (Ad) { 13809566063dSJacob Faibussowitsch PetscCall(PetscObjectBaseTypeCompare((PetscObject)Ad,MATSEQBAIJ,&isbaij)); 13819566063dSJacob Faibussowitsch if (!isbaij) PetscCall(PetscObjectBaseTypeCompare((PetscObject)Ad,MATSEQSBAIJ,&isbaij)); 1382b1a0a8a3SJed Brown } 13838f3b4b4dSDmitry Karpeev if (Ad && nloc > 1) { 1384b1a0a8a3SJed Brown PetscBool match,done; 1385b1a0a8a3SJed Brown /* Try to setup a good matrix partitioning if available */ 13869566063dSJacob Faibussowitsch PetscCall(MatPartitioningCreate(PETSC_COMM_SELF,&mpart)); 13879566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix)); 13889566063dSJacob Faibussowitsch PetscCall(MatPartitioningSetFromOptions(mpart)); 13899566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGCURRENT,&match)); 1390b1a0a8a3SJed Brown if (!match) { 13919566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGSQUARE,&match)); 1392b1a0a8a3SJed Brown } 1393b1a0a8a3SJed Brown if (!match) { /* assume a "good" partitioner is available */ 13941a83f524SJed Brown PetscInt na; 13951a83f524SJed Brown const PetscInt *ia,*ja; 13969566063dSJacob Faibussowitsch PetscCall(MatGetRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done)); 1397b1a0a8a3SJed Brown if (done) { 1398b1a0a8a3SJed Brown /* Build adjacency matrix by hand. Unfortunately a call to 1399b1a0a8a3SJed Brown MatConvert(Ad,MATMPIADJ,MAT_INITIAL_MATRIX,&adj) will 1400b1a0a8a3SJed Brown remove the block-aij structure and we cannot expect 1401b1a0a8a3SJed Brown MatPartitioning to split vertices as we need */ 14020a545947SLisandro Dalcin PetscInt i,j,len,nnz,cnt,*iia=NULL,*jja=NULL; 14031a83f524SJed Brown const PetscInt *row; 1404b1a0a8a3SJed Brown nnz = 0; 1405b1a0a8a3SJed Brown for (i=0; i<na; i++) { /* count number of nonzeros */ 1406b1a0a8a3SJed Brown len = ia[i+1] - ia[i]; 1407b1a0a8a3SJed Brown row = ja + ia[i]; 1408b1a0a8a3SJed Brown for (j=0; j<len; j++) { 1409b1a0a8a3SJed Brown if (row[j] == i) { /* don't count diagonal */ 1410b1a0a8a3SJed Brown len--; break; 1411b1a0a8a3SJed Brown } 1412b1a0a8a3SJed Brown } 1413b1a0a8a3SJed Brown nnz += len; 1414b1a0a8a3SJed Brown } 14159566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(na+1,&iia)); 14169566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nnz,&jja)); 1417b1a0a8a3SJed Brown nnz = 0; 1418b1a0a8a3SJed Brown iia[0] = 0; 1419b1a0a8a3SJed Brown for (i=0; i<na; i++) { /* fill adjacency */ 1420b1a0a8a3SJed Brown cnt = 0; 1421b1a0a8a3SJed Brown len = ia[i+1] - ia[i]; 1422b1a0a8a3SJed Brown row = ja + ia[i]; 1423b1a0a8a3SJed Brown for (j=0; j<len; j++) { 14242fa5cd67SKarl Rupp if (row[j] != i) jja[nnz+cnt++] = row[j]; /* if not diagonal */ 1425b1a0a8a3SJed Brown } 1426b1a0a8a3SJed Brown nnz += cnt; 1427b1a0a8a3SJed Brown iia[i+1] = nnz; 1428b1a0a8a3SJed Brown } 1429b1a0a8a3SJed Brown /* Partitioning of the adjacency matrix */ 14309566063dSJacob Faibussowitsch PetscCall(MatCreateMPIAdj(PETSC_COMM_SELF,na,na,iia,jja,NULL,&adj)); 14319566063dSJacob Faibussowitsch PetscCall(MatPartitioningSetAdjacency(mpart,adj)); 14329566063dSJacob Faibussowitsch PetscCall(MatPartitioningSetNParts(mpart,nloc)); 14339566063dSJacob Faibussowitsch PetscCall(MatPartitioningApply(mpart,&ispart)); 14349566063dSJacob Faibussowitsch PetscCall(ISPartitioningToNumbering(ispart,&isnumb)); 14359566063dSJacob Faibussowitsch PetscCall(MatDestroy(&adj)); 1436b1a0a8a3SJed Brown foundpart = PETSC_TRUE; 1437b1a0a8a3SJed Brown } 14389566063dSJacob Faibussowitsch PetscCall(MatRestoreRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done)); 1439b1a0a8a3SJed Brown } 14409566063dSJacob Faibussowitsch PetscCall(MatPartitioningDestroy(&mpart)); 1441b1a0a8a3SJed Brown } 14429566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nloc,&is)); 1443b1a0a8a3SJed Brown if (!foundpart) { 1444b1a0a8a3SJed Brown 1445b1a0a8a3SJed Brown /* Partitioning by contiguous chunks of rows */ 1446b1a0a8a3SJed Brown 1447b1a0a8a3SJed Brown PetscInt mbs = (rend-rstart)/bs; 1448b1a0a8a3SJed Brown PetscInt start = rstart; 14498f3b4b4dSDmitry Karpeev for (i=0; i<nloc; i++) { 14508f3b4b4dSDmitry Karpeev PetscInt count = (mbs/nloc + ((mbs % nloc) > i)) * bs; 14519566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF,count,start,1,&is[i])); 1452b1a0a8a3SJed Brown start += count; 1453b1a0a8a3SJed Brown } 1454b1a0a8a3SJed Brown 1455b1a0a8a3SJed Brown } else { 1456b1a0a8a3SJed Brown 1457b1a0a8a3SJed Brown /* Partitioning by adjacency of diagonal block */ 1458b1a0a8a3SJed Brown 1459b1a0a8a3SJed Brown const PetscInt *numbering; 1460b1a0a8a3SJed Brown PetscInt *count,nidx,*indices,*newidx,start=0; 1461b1a0a8a3SJed Brown /* Get node count in each partition */ 14629566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nloc,&count)); 14639566063dSJacob Faibussowitsch PetscCall(ISPartitioningCount(ispart,nloc,count)); 1464b1a0a8a3SJed Brown if (isbaij && bs > 1) { /* adjust for the block-aij case */ 14658f3b4b4dSDmitry Karpeev for (i=0; i<nloc; i++) count[i] *= bs; 1466b1a0a8a3SJed Brown } 1467b1a0a8a3SJed Brown /* Build indices from node numbering */ 14689566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(isnumb,&nidx)); 14699566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nidx,&indices)); 1470b1a0a8a3SJed Brown for (i=0; i<nidx; i++) indices[i] = i; /* needs to be initialized */ 14719566063dSJacob Faibussowitsch PetscCall(ISGetIndices(isnumb,&numbering)); 14729566063dSJacob Faibussowitsch PetscCall(PetscSortIntWithPermutation(nidx,numbering,indices)); 14739566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(isnumb,&numbering)); 1474b1a0a8a3SJed Brown if (isbaij && bs > 1) { /* adjust for the block-aij case */ 14759566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nidx*bs,&newidx)); 14762fa5cd67SKarl Rupp for (i=0; i<nidx; i++) { 14772fa5cd67SKarl Rupp for (j=0; j<bs; j++) newidx[i*bs+j] = indices[i]*bs + j; 14782fa5cd67SKarl Rupp } 14799566063dSJacob Faibussowitsch PetscCall(PetscFree(indices)); 1480b1a0a8a3SJed Brown nidx *= bs; 1481b1a0a8a3SJed Brown indices = newidx; 1482b1a0a8a3SJed Brown } 1483b1a0a8a3SJed Brown /* Shift to get global indices */ 1484b1a0a8a3SJed Brown for (i=0; i<nidx; i++) indices[i] += rstart; 1485b1a0a8a3SJed Brown 1486b1a0a8a3SJed Brown /* Build the index sets for each block */ 14878f3b4b4dSDmitry Karpeev for (i=0; i<nloc; i++) { 14889566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF,count[i],&indices[start],PETSC_COPY_VALUES,&is[i])); 14899566063dSJacob Faibussowitsch PetscCall(ISSort(is[i])); 1490b1a0a8a3SJed Brown start += count[i]; 1491b1a0a8a3SJed Brown } 1492b1a0a8a3SJed Brown 14939566063dSJacob Faibussowitsch PetscCall(PetscFree(count)); 14949566063dSJacob Faibussowitsch PetscCall(PetscFree(indices)); 14959566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isnumb)); 14969566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ispart)); 1497b1a0a8a3SJed Brown } 14986a4f0f73SDmitry Karpeev *iis = is; 14998f3b4b4dSDmitry Karpeev PetscFunctionReturn(0); 15008f3b4b4dSDmitry Karpeev } 15018f3b4b4dSDmitry Karpeev 1502b541e6a4SDmitry Karpeev PETSC_INTERN PetscErrorCode PCGASMCreateStraddlingSubdomains(Mat A,PetscInt N,PetscInt *n,IS *iis[]) 15038f3b4b4dSDmitry Karpeev { 15048f3b4b4dSDmitry Karpeev PetscFunctionBegin; 15059566063dSJacob Faibussowitsch PetscCall(MatSubdomainsCreateCoalesce(A,N,n,iis)); 15068f3b4b4dSDmitry Karpeev PetscFunctionReturn(0); 15078f3b4b4dSDmitry Karpeev } 15088f3b4b4dSDmitry Karpeev 15098f3b4b4dSDmitry Karpeev /*@C 15108f3b4b4dSDmitry Karpeev PCGASMCreateSubdomains - Creates n index sets defining n nonoverlapping subdomains for the additive 15118f3b4b4dSDmitry Karpeev Schwarz preconditioner for a any problem based on its matrix. 15128f3b4b4dSDmitry Karpeev 15138f3b4b4dSDmitry Karpeev Collective 15148f3b4b4dSDmitry Karpeev 15158f3b4b4dSDmitry Karpeev Input Parameters: 15168f3b4b4dSDmitry Karpeev + A - The global matrix operator 15178f3b4b4dSDmitry Karpeev - N - the number of global subdomains requested 15188f3b4b4dSDmitry Karpeev 15198f3b4b4dSDmitry Karpeev Output Parameters: 15208f3b4b4dSDmitry Karpeev + n - the number of subdomains created on this processor 15218f3b4b4dSDmitry Karpeev - iis - the array of index sets defining the local inner subdomains (on which the correction is applied) 15228f3b4b4dSDmitry Karpeev 15238f3b4b4dSDmitry Karpeev Level: advanced 15248f3b4b4dSDmitry Karpeev 15258f3b4b4dSDmitry Karpeev Note: When N >= A's communicator size, each subdomain is local -- contained within a single processor. 15268f3b4b4dSDmitry Karpeev When N < size, the subdomains are 'straddling' (processor boundaries) and are no longer local. 15278f3b4b4dSDmitry Karpeev The resulting subdomains can be use in PCGASMSetSubdomains(pc,n,iss,NULL). The overlapping 15288f3b4b4dSDmitry Karpeev outer subdomains will be automatically generated from these according to the requested amount of 15298f3b4b4dSDmitry Karpeev overlap; this is currently supported only with local subdomains. 15308f3b4b4dSDmitry Karpeev 1531db781477SPatrick Sanan .seealso: `PCGASMSetSubdomains()`, `PCGASMDestroySubdomains()` 15328f3b4b4dSDmitry Karpeev @*/ 1533b541e6a4SDmitry Karpeev PetscErrorCode PCGASMCreateSubdomains(Mat A,PetscInt N,PetscInt *n,IS *iis[]) 15348f3b4b4dSDmitry Karpeev { 15358f3b4b4dSDmitry Karpeev PetscMPIInt size; 15368f3b4b4dSDmitry Karpeev 15378f3b4b4dSDmitry Karpeev PetscFunctionBegin; 15388f3b4b4dSDmitry Karpeev PetscValidHeaderSpecific(A,MAT_CLASSID,1); 15398f3b4b4dSDmitry Karpeev PetscValidPointer(iis,4); 15408f3b4b4dSDmitry Karpeev 154163a3b9bcSJacob Faibussowitsch PetscCheck(N >= 1,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Number of subdomains must be > 0, N = %" PetscInt_FMT,N); 15429566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A),&size)); 15438f3b4b4dSDmitry Karpeev if (N >= size) { 15448f3b4b4dSDmitry Karpeev *n = N/size + (N%size); 15459566063dSJacob Faibussowitsch PetscCall(PCGASMCreateLocalSubdomains(A,*n,iis)); 15466a4f0f73SDmitry Karpeev } else { 15479566063dSJacob Faibussowitsch PetscCall(PCGASMCreateStraddlingSubdomains(A,N,n,iis)); 15486a4f0f73SDmitry Karpeev } 1549b1a0a8a3SJed Brown PetscFunctionReturn(0); 1550b1a0a8a3SJed Brown } 1551b1a0a8a3SJed Brown 1552b1a0a8a3SJed Brown /*@C 1553f746d493SDmitry Karpeev PCGASMDestroySubdomains - Destroys the index sets created with 15548f3b4b4dSDmitry Karpeev PCGASMCreateSubdomains() or PCGASMCreateSubdomains2D. Should be 155506b43e7eSDmitry Karpeev called after setting subdomains with PCGASMSetSubdomains(). 1556b1a0a8a3SJed Brown 1557b1a0a8a3SJed Brown Collective 1558b1a0a8a3SJed Brown 1559b1a0a8a3SJed Brown Input Parameters: 1560b1a0a8a3SJed Brown + n - the number of index sets 15616a4f0f73SDmitry Karpeev . iis - the array of inner subdomains, 15620298fd71SBarry Smith - ois - the array of outer subdomains, can be NULL 1563b1a0a8a3SJed Brown 15646a4f0f73SDmitry Karpeev Level: intermediate 15656a4f0f73SDmitry Karpeev 156695452b02SPatrick Sanan Notes: 156795452b02SPatrick Sanan this is merely a convenience subroutine that walks each list, 15682c112581SDmitry Karpeev destroys each IS on the list, and then frees the list. At the end the 15692c112581SDmitry Karpeev list pointers are set to NULL. 1570b1a0a8a3SJed Brown 1571db781477SPatrick Sanan .seealso: `PCGASMCreateSubdomains()`, `PCGASMSetSubdomains()` 1572b1a0a8a3SJed Brown @*/ 15732c112581SDmitry Karpeev PetscErrorCode PCGASMDestroySubdomains(PetscInt n,IS **iis,IS **ois) 1574b1a0a8a3SJed Brown { 1575b1a0a8a3SJed Brown PetscInt i; 15765fd66863SKarl Rupp 1577b1a0a8a3SJed Brown PetscFunctionBegin; 1578a35b7b57SDmitry Karpeev if (n <= 0) PetscFunctionReturn(0); 15796a4f0f73SDmitry Karpeev if (ois) { 15802c112581SDmitry Karpeev PetscValidPointer(ois,3); 15812c112581SDmitry Karpeev if (*ois) { 15822c112581SDmitry Karpeev PetscValidPointer(*ois,3); 1583a35b7b57SDmitry Karpeev for (i=0; i<n; i++) { 15849566063dSJacob Faibussowitsch PetscCall(ISDestroy(&(*ois)[i])); 1585a35b7b57SDmitry Karpeev } 15869566063dSJacob Faibussowitsch PetscCall(PetscFree((*ois))); 15872c112581SDmitry Karpeev } 1588b1a0a8a3SJed Brown } 1589b9d0fdaaSFande Kong if (iis) { 1590b9d0fdaaSFande Kong PetscValidPointer(iis,2); 1591b9d0fdaaSFande Kong if (*iis) { 1592b9d0fdaaSFande Kong PetscValidPointer(*iis,2); 1593b9d0fdaaSFande Kong for (i=0; i<n; i++) { 15949566063dSJacob Faibussowitsch PetscCall(ISDestroy(&(*iis)[i])); 1595b9d0fdaaSFande Kong } 15969566063dSJacob Faibussowitsch PetscCall(PetscFree((*iis))); 1597b9d0fdaaSFande Kong } 1598b9d0fdaaSFande Kong } 1599b1a0a8a3SJed Brown PetscFunctionReturn(0); 1600b1a0a8a3SJed Brown } 1601b1a0a8a3SJed Brown 1602af538404SDmitry Karpeev #define PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,xleft_loc,ylow_loc,xright_loc,yhigh_loc,n) \ 1603af538404SDmitry Karpeev { \ 1604af538404SDmitry Karpeev PetscInt first_row = first/M, last_row = last/M+1; \ 1605af538404SDmitry Karpeev /* \ 1606af538404SDmitry Karpeev Compute ylow_loc and yhigh_loc so that (ylow_loc,xleft) and (yhigh_loc,xright) are the corners \ 1607af538404SDmitry Karpeev of the bounding box of the intersection of the subdomain with the local ownership range (local \ 1608af538404SDmitry Karpeev subdomain). \ 1609af538404SDmitry Karpeev Also compute xleft_loc and xright_loc as the lower and upper bounds on the first and last rows \ 1610af538404SDmitry Karpeev of the intersection. \ 1611af538404SDmitry Karpeev */ \ 1612af538404SDmitry Karpeev /* ylow_loc is the grid row containing the first element of the local sumbdomain */ \ 1613eec7e3faSDmitry Karpeev *ylow_loc = PetscMax(first_row,ylow); \ 1614af538404SDmitry Karpeev /* xleft_loc is the offset of first element of the local subdomain within its grid row (might actually be outside the local subdomain) */ \ 1615eec7e3faSDmitry Karpeev *xleft_loc = *ylow_loc==first_row ? PetscMax(first%M,xleft) : xleft; \ 1616af538404SDmitry Karpeev /* yhigh_loc is the grid row above the last local subdomain element */ \ 1617eec7e3faSDmitry Karpeev *yhigh_loc = PetscMin(last_row,yhigh); \ 1618af538404SDmitry Karpeev /* xright is the offset of the end of the local subdomain within its grid row (might actually be outside the local subdomain) */ \ 1619eec7e3faSDmitry Karpeev *xright_loc = *yhigh_loc==last_row ? PetscMin(xright,last%M) : xright; \ 1620af538404SDmitry Karpeev /* Now compute the size of the local subdomain n. */ \ 1621c3518bceSDmitry Karpeev *n = 0; \ 1622eec7e3faSDmitry Karpeev if (*ylow_loc < *yhigh_loc) { \ 1623af538404SDmitry Karpeev PetscInt width = xright-xleft; \ 1624eec7e3faSDmitry Karpeev *n += width*(*yhigh_loc-*ylow_loc-1); \ 1625eec7e3faSDmitry Karpeev *n += PetscMin(PetscMax(*xright_loc-xleft,0),width); \ 1626eec7e3faSDmitry Karpeev *n -= PetscMin(PetscMax(*xleft_loc-xleft,0), width); \ 1627af538404SDmitry Karpeev } \ 1628af538404SDmitry Karpeev } 1629af538404SDmitry Karpeev 1630b1a0a8a3SJed Brown /*@ 1631f746d493SDmitry Karpeev PCGASMCreateSubdomains2D - Creates the index sets for the overlapping Schwarz 1632b1a0a8a3SJed Brown preconditioner for a two-dimensional problem on a regular grid. 1633b1a0a8a3SJed Brown 1634af538404SDmitry Karpeev Collective 1635b1a0a8a3SJed Brown 1636b1a0a8a3SJed Brown Input Parameters: 16376b867d5aSJose E. Roman + pc - the preconditioner context 16386b867d5aSJose E. Roman . M - the global number of grid points in the x direction 16396b867d5aSJose E. Roman . N - the global number of grid points in the y direction 16406b867d5aSJose E. Roman . Mdomains - the global number of subdomains in the x direction 16416b867d5aSJose E. Roman . Ndomains - the global number of subdomains in the y direction 1642b1a0a8a3SJed Brown . dof - degrees of freedom per node 1643b1a0a8a3SJed Brown - overlap - overlap in mesh lines 1644b1a0a8a3SJed Brown 1645b1a0a8a3SJed Brown Output Parameters: 1646af538404SDmitry Karpeev + Nsub - the number of local subdomains created 16476a4f0f73SDmitry Karpeev . iis - array of index sets defining inner (nonoverlapping) subdomains 16486a4f0f73SDmitry Karpeev - ois - array of index sets defining outer (overlapping, if overlap > 0) subdomains 1649b1a0a8a3SJed Brown 1650b1a0a8a3SJed Brown Level: advanced 1651b1a0a8a3SJed Brown 1652db781477SPatrick Sanan .seealso: `PCGASMSetSubdomains()`, `PCGASMGetSubKSP()`, `PCGASMSetOverlap()` 1653b1a0a8a3SJed Brown @*/ 16546a4f0f73SDmitry Karpeev PetscErrorCode PCGASMCreateSubdomains2D(PC pc,PetscInt M,PetscInt N,PetscInt Mdomains,PetscInt Ndomains,PetscInt dof,PetscInt overlap,PetscInt *nsub,IS **iis,IS **ois) 1655b1a0a8a3SJed Brown { 16562bbb417fSDmitry Karpeev PetscMPIInt size, rank; 16572bbb417fSDmitry Karpeev PetscInt i, j; 16582bbb417fSDmitry Karpeev PetscInt maxheight, maxwidth; 16592bbb417fSDmitry Karpeev PetscInt xstart, xleft, xright, xleft_loc, xright_loc; 16602bbb417fSDmitry Karpeev PetscInt ystart, ylow, yhigh, ylow_loc, yhigh_loc; 1661eec7e3faSDmitry Karpeev PetscInt x[2][2], y[2][2], n[2]; 16622bbb417fSDmitry Karpeev PetscInt first, last; 16632bbb417fSDmitry Karpeev PetscInt nidx, *idx; 16642bbb417fSDmitry Karpeev PetscInt ii,jj,s,q,d; 1665af538404SDmitry Karpeev PetscInt k,kk; 16662bbb417fSDmitry Karpeev PetscMPIInt color; 16672bbb417fSDmitry Karpeev MPI_Comm comm, subcomm; 16680a545947SLisandro Dalcin IS **xis = NULL, **is = ois, **is_local = iis; 1669d34fcf5fSBarry Smith 1670b1a0a8a3SJed Brown PetscFunctionBegin; 16719566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)pc, &comm)); 16729566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 16739566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 16749566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(pc->pmat, &first, &last)); 16752472a847SBarry Smith PetscCheck((first%dof) == 0 && (last%dof) == 0,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Matrix row partitioning unsuitable for domain decomposition: local row range (%" PetscInt_FMT ",%" PetscInt_FMT ") " 167663a3b9bcSJacob Faibussowitsch "does not respect the number of degrees of freedom per grid point %" PetscInt_FMT, first, last, dof); 1677d34fcf5fSBarry Smith 1678af538404SDmitry Karpeev /* Determine the number of domains with nonzero intersections with the local ownership range. */ 16792bbb417fSDmitry Karpeev s = 0; 1680af538404SDmitry Karpeev ystart = 0; 1681af538404SDmitry Karpeev for (j=0; j<Ndomains; ++j) { 1682af538404SDmitry Karpeev maxheight = N/Ndomains + ((N % Ndomains) > j); /* Maximal height of subdomain */ 168363a3b9bcSJacob Faibussowitsch PetscCheck(maxheight >= 2,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many %" PetscInt_FMT " subdomains in the vertical directon for mesh height %" PetscInt_FMT, Ndomains, N); 1684eec7e3faSDmitry Karpeev /* Vertical domain limits with an overlap. */ 1685eec7e3faSDmitry Karpeev ylow = PetscMax(ystart - overlap,0); 1686eec7e3faSDmitry Karpeev yhigh = PetscMin(ystart + maxheight + overlap,N); 1687b1a0a8a3SJed Brown xstart = 0; 1688af538404SDmitry Karpeev for (i=0; i<Mdomains; ++i) { 1689af538404SDmitry Karpeev maxwidth = M/Mdomains + ((M % Mdomains) > i); /* Maximal width of subdomain */ 169063a3b9bcSJacob Faibussowitsch PetscCheck(maxwidth >= 2,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many %" PetscInt_FMT " subdomains in the horizontal direction for mesh width %" PetscInt_FMT, Mdomains, M); 1691eec7e3faSDmitry Karpeev /* Horizontal domain limits with an overlap. */ 1692eec7e3faSDmitry Karpeev xleft = PetscMax(xstart - overlap,0); 1693eec7e3faSDmitry Karpeev xright = PetscMin(xstart + maxwidth + overlap,M); 1694af538404SDmitry Karpeev /* 1695af538404SDmitry Karpeev Determine whether this subdomain intersects this processor's ownership range of pc->pmat. 1696af538404SDmitry Karpeev */ 1697c3518bceSDmitry Karpeev PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,(&xleft_loc),(&ylow_loc),(&xright_loc),(&yhigh_loc),(&nidx)); 16982fa5cd67SKarl Rupp if (nidx) ++s; 1699af538404SDmitry Karpeev xstart += maxwidth; 1700af538404SDmitry Karpeev } /* for (i = 0; i < Mdomains; ++i) */ 1701af538404SDmitry Karpeev ystart += maxheight; 1702af538404SDmitry Karpeev } /* for (j = 0; j < Ndomains; ++j) */ 17032fa5cd67SKarl Rupp 1704af538404SDmitry Karpeev /* Now we can allocate the necessary number of ISs. */ 1705af538404SDmitry Karpeev *nsub = s; 17069566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(*nsub,is)); 17079566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(*nsub,is_local)); 1708af538404SDmitry Karpeev s = 0; 1709af538404SDmitry Karpeev ystart = 0; 1710af538404SDmitry Karpeev for (j=0; j<Ndomains; ++j) { 1711af538404SDmitry Karpeev maxheight = N/Ndomains + ((N % Ndomains) > j); /* Maximal height of subdomain */ 171263a3b9bcSJacob Faibussowitsch PetscCheck(maxheight >= 2,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many %" PetscInt_FMT " subdomains in the vertical directon for mesh height %" PetscInt_FMT, Ndomains, N); 1713af538404SDmitry Karpeev /* Vertical domain limits with an overlap. */ 1714eec7e3faSDmitry Karpeev y[0][0] = PetscMax(ystart - overlap,0); 1715eec7e3faSDmitry Karpeev y[0][1] = PetscMin(ystart + maxheight + overlap,N); 1716eec7e3faSDmitry Karpeev /* Vertical domain limits without an overlap. */ 1717eec7e3faSDmitry Karpeev y[1][0] = ystart; 1718eec7e3faSDmitry Karpeev y[1][1] = PetscMin(ystart + maxheight,N); 1719eec7e3faSDmitry Karpeev xstart = 0; 1720af538404SDmitry Karpeev for (i=0; i<Mdomains; ++i) { 1721af538404SDmitry Karpeev maxwidth = M/Mdomains + ((M % Mdomains) > i); /* Maximal width of subdomain */ 172263a3b9bcSJacob Faibussowitsch PetscCheck(maxwidth >= 2,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many %" PetscInt_FMT " subdomains in the horizontal direction for mesh width %" PetscInt_FMT, Mdomains, M); 1723af538404SDmitry Karpeev /* Horizontal domain limits with an overlap. */ 1724eec7e3faSDmitry Karpeev x[0][0] = PetscMax(xstart - overlap,0); 1725eec7e3faSDmitry Karpeev x[0][1] = PetscMin(xstart + maxwidth + overlap,M); 1726eec7e3faSDmitry Karpeev /* Horizontal domain limits without an overlap. */ 1727eec7e3faSDmitry Karpeev x[1][0] = xstart; 1728eec7e3faSDmitry Karpeev x[1][1] = PetscMin(xstart+maxwidth,M); 17292bbb417fSDmitry Karpeev /* 17302bbb417fSDmitry Karpeev Determine whether this domain intersects this processor's ownership range of pc->pmat. 17312bbb417fSDmitry Karpeev Do this twice: first for the domains with overlaps, and once without. 17322bbb417fSDmitry Karpeev During the first pass create the subcommunicators, and use them on the second pass as well. 17332bbb417fSDmitry Karpeev */ 17342bbb417fSDmitry Karpeev for (q = 0; q < 2; ++q) { 1735cc96b2e5SBarry Smith PetscBool split = PETSC_FALSE; 17362bbb417fSDmitry Karpeev /* 17372bbb417fSDmitry Karpeev domain limits, (xleft, xright) and (ylow, yheigh) are adjusted 17382bbb417fSDmitry Karpeev according to whether the domain with an overlap or without is considered. 17392bbb417fSDmitry Karpeev */ 17402bbb417fSDmitry Karpeev xleft = x[q][0]; xright = x[q][1]; 17412bbb417fSDmitry Karpeev ylow = y[q][0]; yhigh = y[q][1]; 1742c3518bceSDmitry Karpeev PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,(&xleft_loc),(&ylow_loc),(&xright_loc),(&yhigh_loc),(&nidx)); 1743af538404SDmitry Karpeev nidx *= dof; 1744eec7e3faSDmitry Karpeev n[q] = nidx; 17452bbb417fSDmitry Karpeev /* 1746eec7e3faSDmitry Karpeev Based on the counted number of indices in the local domain *with an overlap*, 1747af538404SDmitry Karpeev construct a subcommunicator of all the processors supporting this domain. 1748eec7e3faSDmitry Karpeev Observe that a domain with an overlap might have nontrivial local support, 1749eec7e3faSDmitry Karpeev while the domain without an overlap might not. Hence, the decision to participate 1750eec7e3faSDmitry Karpeev in the subcommunicator must be based on the domain with an overlap. 17512bbb417fSDmitry Karpeev */ 17522bbb417fSDmitry Karpeev if (q == 0) { 17532fa5cd67SKarl Rupp if (nidx) color = 1; 17542fa5cd67SKarl Rupp else color = MPI_UNDEFINED; 17559566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_split(comm, color, rank, &subcomm)); 1756cc96b2e5SBarry Smith split = PETSC_TRUE; 17572bbb417fSDmitry Karpeev } 1758af538404SDmitry Karpeev /* 1759eec7e3faSDmitry Karpeev Proceed only if the number of local indices *with an overlap* is nonzero. 1760af538404SDmitry Karpeev */ 1761eec7e3faSDmitry Karpeev if (n[0]) { 17622fa5cd67SKarl Rupp if (q == 0) xis = is; 1763af538404SDmitry Karpeev if (q == 1) { 1764af538404SDmitry Karpeev /* 1765eec7e3faSDmitry Karpeev The IS for the no-overlap subdomain shares a communicator with the overlapping domain. 1766af538404SDmitry Karpeev Moreover, if the overlap is zero, the two ISs are identical. 1767af538404SDmitry Karpeev */ 1768b1a0a8a3SJed Brown if (overlap == 0) { 1769eec7e3faSDmitry Karpeev (*is_local)[s] = (*is)[s]; 17709566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)(*is)[s])); 17712bbb417fSDmitry Karpeev continue; 1772d34fcf5fSBarry Smith } else { 17736a4f0f73SDmitry Karpeev xis = is_local; 1774eec7e3faSDmitry Karpeev subcomm = ((PetscObject)(*is)[s])->comm; 17752bbb417fSDmitry Karpeev } 1776af538404SDmitry Karpeev } /* if (q == 1) */ 17770298fd71SBarry Smith idx = NULL; 17789566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nidx,&idx)); 1779eec7e3faSDmitry Karpeev if (nidx) { 17802bbb417fSDmitry Karpeev k = 0; 17812bbb417fSDmitry Karpeev for (jj=ylow_loc; jj<yhigh_loc; ++jj) { 1782af538404SDmitry Karpeev PetscInt x0 = (jj==ylow_loc) ? xleft_loc : xleft; 1783af538404SDmitry Karpeev PetscInt x1 = (jj==yhigh_loc-1) ? xright_loc : xright; 1784af538404SDmitry Karpeev kk = dof*(M*jj + x0); 17852bbb417fSDmitry Karpeev for (ii=x0; ii<x1; ++ii) { 17862bbb417fSDmitry Karpeev for (d = 0; d < dof; ++d) { 17872bbb417fSDmitry Karpeev idx[k++] = kk++; 1788b1a0a8a3SJed Brown } 1789b1a0a8a3SJed Brown } 1790b1a0a8a3SJed Brown } 1791eec7e3faSDmitry Karpeev } 17929566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(subcomm,nidx,idx,PETSC_OWN_POINTER,(*xis)+s)); 1793cc96b2e5SBarry Smith if (split) { 17949566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free(&subcomm)); 1795cc96b2e5SBarry Smith } 1796eec7e3faSDmitry Karpeev }/* if (n[0]) */ 17972bbb417fSDmitry Karpeev }/* for (q = 0; q < 2; ++q) */ 17982fa5cd67SKarl Rupp if (n[0]) ++s; 1799af538404SDmitry Karpeev xstart += maxwidth; 1800af538404SDmitry Karpeev } /* for (i = 0; i < Mdomains; ++i) */ 18012bbb417fSDmitry Karpeev ystart += maxheight; 1802af538404SDmitry Karpeev } /* for (j = 0; j < Ndomains; ++j) */ 1803b1a0a8a3SJed Brown PetscFunctionReturn(0); 1804b1a0a8a3SJed Brown } 1805b1a0a8a3SJed Brown 1806b1a0a8a3SJed Brown /*@C 180706b43e7eSDmitry Karpeev PCGASMGetSubdomains - Gets the subdomains supported on this processor 180806b43e7eSDmitry Karpeev for the additive Schwarz preconditioner. 1809b1a0a8a3SJed Brown 1810b1a0a8a3SJed Brown Not Collective 1811b1a0a8a3SJed Brown 1812b1a0a8a3SJed Brown Input Parameter: 1813b1a0a8a3SJed Brown . pc - the preconditioner context 1814b1a0a8a3SJed Brown 1815b1a0a8a3SJed Brown Output Parameters: 1816b1a0a8a3SJed Brown + n - the number of subdomains for this processor (default value = 1) 18170298fd71SBarry Smith . iis - the index sets that define the inner subdomains (without overlap) supported on this processor (can be NULL) 18180298fd71SBarry Smith - ois - the index sets that define the outer subdomains (with overlap) supported on this processor (can be NULL) 1819b1a0a8a3SJed Brown 1820b1a0a8a3SJed Brown Notes: 18216a4f0f73SDmitry Karpeev The user is responsible for destroying the ISs and freeing the returned arrays. 1822b1a0a8a3SJed Brown The IS numbering is in the parallel, global numbering of the vector. 1823b1a0a8a3SJed Brown 1824b1a0a8a3SJed Brown Level: advanced 1825b1a0a8a3SJed Brown 1826db781477SPatrick Sanan .seealso: `PCGASMSetOverlap()`, `PCGASMGetSubKSP()`, `PCGASMCreateSubdomains2D()`, 1827db781477SPatrick Sanan `PCGASMSetSubdomains()`, `PCGASMGetSubmatrices()` 1828b1a0a8a3SJed Brown @*/ 18296a4f0f73SDmitry Karpeev PetscErrorCode PCGASMGetSubdomains(PC pc,PetscInt *n,IS *iis[],IS *ois[]) 1830b1a0a8a3SJed Brown { 1831f746d493SDmitry Karpeev PC_GASM *osm; 1832b1a0a8a3SJed Brown PetscBool match; 18336a4f0f73SDmitry Karpeev PetscInt i; 18345fd66863SKarl Rupp 1835b1a0a8a3SJed Brown PetscFunctionBegin; 1836b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 18379566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match)); 183828b400f6SJacob Faibussowitsch PetscCheck(match,PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Incorrect object type: expected %s, got %s instead", PCGASM, ((PetscObject)pc)->type_name); 1839f746d493SDmitry Karpeev osm = (PC_GASM*)pc->data; 1840ab3e923fSDmitry Karpeev if (n) *n = osm->n; 1841*1baa6e33SBarry Smith if (iis) PetscCall(PetscMalloc1(osm->n, iis)); 1842*1baa6e33SBarry Smith if (ois) PetscCall(PetscMalloc1(osm->n, ois)); 18436a4f0f73SDmitry Karpeev if (iis || ois) { 18446a4f0f73SDmitry Karpeev for (i = 0; i < osm->n; ++i) { 18456a4f0f73SDmitry Karpeev if (iis) (*iis)[i] = osm->iis[i]; 18466a4f0f73SDmitry Karpeev if (ois) (*ois)[i] = osm->ois[i]; 18476a4f0f73SDmitry Karpeev } 1848b1a0a8a3SJed Brown } 1849b1a0a8a3SJed Brown PetscFunctionReturn(0); 1850b1a0a8a3SJed Brown } 1851b1a0a8a3SJed Brown 1852b1a0a8a3SJed Brown /*@C 185306b43e7eSDmitry Karpeev PCGASMGetSubmatrices - Gets the local submatrices (for this processor 1854b1a0a8a3SJed Brown only) for the additive Schwarz preconditioner. 1855b1a0a8a3SJed Brown 1856b1a0a8a3SJed Brown Not Collective 1857b1a0a8a3SJed Brown 1858b1a0a8a3SJed Brown Input Parameter: 1859b1a0a8a3SJed Brown . pc - the preconditioner context 1860b1a0a8a3SJed Brown 1861b1a0a8a3SJed Brown Output Parameters: 1862b1a0a8a3SJed Brown + n - the number of matrices for this processor (default value = 1) 1863b1a0a8a3SJed Brown - mat - the matrices 1864b1a0a8a3SJed Brown 186595452b02SPatrick Sanan Notes: 186695452b02SPatrick Sanan matrices returned by this routine have the same communicators as the index sets (IS) 18678f3b4b4dSDmitry Karpeev used to define subdomains in PCGASMSetSubdomains() 1868b1a0a8a3SJed Brown Level: advanced 1869b1a0a8a3SJed Brown 1870db781477SPatrick Sanan .seealso: `PCGASMSetOverlap()`, `PCGASMGetSubKSP()`, 1871db781477SPatrick Sanan `PCGASMCreateSubdomains2D()`, `PCGASMSetSubdomains()`, `PCGASMGetSubdomains()` 1872b1a0a8a3SJed Brown @*/ 187306b43e7eSDmitry Karpeev PetscErrorCode PCGASMGetSubmatrices(PC pc,PetscInt *n,Mat *mat[]) 1874b1a0a8a3SJed Brown { 1875f746d493SDmitry Karpeev PC_GASM *osm; 1876b1a0a8a3SJed Brown PetscBool match; 1877b1a0a8a3SJed Brown 1878b1a0a8a3SJed Brown PetscFunctionBegin; 1879b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1880b1a0a8a3SJed Brown PetscValidIntPointer(n,2); 1881b1a0a8a3SJed Brown if (mat) PetscValidPointer(mat,3); 188228b400f6SJacob Faibussowitsch PetscCheck(pc->setupcalled,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call after KSPSetUp() or PCSetUp()."); 18839566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match)); 188428b400f6SJacob Faibussowitsch PetscCheck(match,PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Expected %s, got %s instead", PCGASM, ((PetscObject)pc)->type_name); 1885f746d493SDmitry Karpeev osm = (PC_GASM*)pc->data; 1886ab3e923fSDmitry Karpeev if (n) *n = osm->n; 1887b1a0a8a3SJed Brown if (mat) *mat = osm->pmat; 1888b1a0a8a3SJed Brown PetscFunctionReturn(0); 1889b1a0a8a3SJed Brown } 1890d709ea83SDmitry Karpeev 1891d709ea83SDmitry Karpeev /*@ 18928f3b4b4dSDmitry Karpeev PCGASMSetUseDMSubdomains - Indicates whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible. 1893d709ea83SDmitry Karpeev Logically Collective 1894d709ea83SDmitry Karpeev 1895d8d19677SJose E. Roman Input Parameters: 1896d709ea83SDmitry Karpeev + pc - the preconditioner 1897d709ea83SDmitry Karpeev - flg - boolean indicating whether to use subdomains defined by the DM 1898d709ea83SDmitry Karpeev 1899d709ea83SDmitry Karpeev Options Database Key: 19008f3b4b4dSDmitry Karpeev . -pc_gasm_dm_subdomains -pc_gasm_overlap -pc_gasm_total_subdomains 1901d709ea83SDmitry Karpeev 1902d709ea83SDmitry Karpeev Level: intermediate 1903d709ea83SDmitry Karpeev 1904d709ea83SDmitry Karpeev Notes: 19058f3b4b4dSDmitry Karpeev PCGASMSetSubdomains(), PCGASMSetTotalSubdomains() or PCGASMSetOverlap() take precedence over PCGASMSetUseDMSubdomains(), 19068f3b4b4dSDmitry Karpeev so setting PCGASMSetSubdomains() with nontrivial subdomain ISs or any of PCGASMSetTotalSubdomains() and PCGASMSetOverlap() 19078f3b4b4dSDmitry Karpeev automatically turns the latter off. 1908d709ea83SDmitry Karpeev 1909db781477SPatrick Sanan .seealso: `PCGASMGetUseDMSubdomains()`, `PCGASMSetSubdomains()`, `PCGASMSetOverlap()` 1910db781477SPatrick Sanan `PCGASMCreateSubdomains2D()` 1911d709ea83SDmitry Karpeev @*/ 19128f3b4b4dSDmitry Karpeev PetscErrorCode PCGASMSetUseDMSubdomains(PC pc,PetscBool flg) 1913d709ea83SDmitry Karpeev { 1914d709ea83SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 1915d709ea83SDmitry Karpeev PetscBool match; 1916d709ea83SDmitry Karpeev 1917d709ea83SDmitry Karpeev PetscFunctionBegin; 1918d709ea83SDmitry Karpeev PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1919d709ea83SDmitry Karpeev PetscValidLogicalCollectiveBool(pc,flg,2); 192028b400f6SJacob Faibussowitsch PetscCheck(!pc->setupcalled,((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Not for a setup PC."); 19219566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match)); 1922d709ea83SDmitry Karpeev if (match) { 19238f3b4b4dSDmitry Karpeev if (!osm->user_subdomains && osm->N == PETSC_DETERMINE && osm->overlap < 0) { 1924d709ea83SDmitry Karpeev osm->dm_subdomains = flg; 1925d709ea83SDmitry Karpeev } 19268f3b4b4dSDmitry Karpeev } 1927d709ea83SDmitry Karpeev PetscFunctionReturn(0); 1928d709ea83SDmitry Karpeev } 1929d709ea83SDmitry Karpeev 1930d709ea83SDmitry Karpeev /*@ 19318f3b4b4dSDmitry Karpeev PCGASMGetUseDMSubdomains - Returns flag indicating whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible. 1932d709ea83SDmitry Karpeev Not Collective 1933d709ea83SDmitry Karpeev 1934d709ea83SDmitry Karpeev Input Parameter: 1935d709ea83SDmitry Karpeev . pc - the preconditioner 1936d709ea83SDmitry Karpeev 1937d709ea83SDmitry Karpeev Output Parameter: 1938d709ea83SDmitry Karpeev . flg - boolean indicating whether to use subdomains defined by the DM 1939d709ea83SDmitry Karpeev 1940d709ea83SDmitry Karpeev Level: intermediate 1941d709ea83SDmitry Karpeev 1942db781477SPatrick Sanan .seealso: `PCGASMSetUseDMSubdomains()`, `PCGASMSetOverlap()` 1943db781477SPatrick Sanan `PCGASMCreateSubdomains2D()` 1944d709ea83SDmitry Karpeev @*/ 19458f3b4b4dSDmitry Karpeev PetscErrorCode PCGASMGetUseDMSubdomains(PC pc,PetscBool* flg) 1946d709ea83SDmitry Karpeev { 1947d709ea83SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 1948d709ea83SDmitry Karpeev PetscBool match; 1949d709ea83SDmitry Karpeev 1950d709ea83SDmitry Karpeev PetscFunctionBegin; 1951d709ea83SDmitry Karpeev PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1952534a8f05SLisandro Dalcin PetscValidBoolPointer(flg,2); 19539566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match)); 1954d709ea83SDmitry Karpeev if (match) { 1955d709ea83SDmitry Karpeev if (flg) *flg = osm->dm_subdomains; 1956d709ea83SDmitry Karpeev } 1957d709ea83SDmitry Karpeev PetscFunctionReturn(0); 1958d709ea83SDmitry Karpeev } 1959