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)); 22906b43e7eSDmitry Karpeev if (view_subdomains) { 2309566063dSJacob Faibussowitsch PetscCall(PCGASMSubdomainView_Private(pc,d,sviewer)); 23106b43e7eSDmitry Karpeev } 2326a4f0f73SDmitry Karpeev if (!pc->setupcalled) { 2339566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(sviewer, " Solver not set up yet: PCSetUp() not yet called\n")); 2342fa5cd67SKarl Rupp } else { 2359566063dSJacob Faibussowitsch PetscCall(KSPView(osm->ksp[d],sviewer)); 2366a4f0f73SDmitry Karpeev } 2379566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(sviewer," - - - - - - - - - - - - - - - - - -\n")); 2389566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(sviewer)); 2399566063dSJacob Faibussowitsch PetscCall(PetscViewerRestoreSubViewer(viewer,((PetscObject)osm->ois[d])->comm, &sviewer)); 24006b43e7eSDmitry Karpeev ++l; 241b1a0a8a3SJed Brown } 24206b43e7eSDmitry Karpeev } 2439566063dSJacob Faibussowitsch PetscCallMPI(MPI_Barrier(PetscObjectComm((PetscObject)pc))); 244b1a0a8a3SJed Brown } 2459566063dSJacob Faibussowitsch PetscCall(PetscFree2(numbering,permutation)); 2469566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 2479566063dSJacob Faibussowitsch PetscCall(PetscViewerFlush(viewer)); 2489530cbd7SBarry Smith /* this line is needed to match the extra PetscViewerASCIIPushSynchronized() in PetscViewerGetSubViewer() */ 2499566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopSynchronized(viewer)); 2509530cbd7SBarry Smith } 251b1a0a8a3SJed Brown PetscFunctionReturn(0); 252b1a0a8a3SJed Brown } 253b1a0a8a3SJed Brown 2548f3b4b4dSDmitry Karpeev PETSC_INTERN PetscErrorCode PCGASMCreateLocalSubdomains(Mat A, PetscInt nloc, IS *iis[]); 255af538404SDmitry Karpeev 256ea91fabdSFande Kong PetscErrorCode PCGASMSetHierarchicalPartitioning(PC pc) 257ea91fabdSFande Kong { 258ea91fabdSFande Kong PC_GASM *osm = (PC_GASM*)pc->data; 259ea91fabdSFande Kong MatPartitioning part; 260ea91fabdSFande Kong MPI_Comm comm; 261ea91fabdSFande Kong PetscMPIInt size; 262ea91fabdSFande Kong PetscInt nlocalsubdomains,fromrows_localsize; 263ea91fabdSFande Kong IS partitioning,fromrows,isn; 264ea91fabdSFande Kong Vec outervec; 265ea91fabdSFande Kong 266ea91fabdSFande Kong PetscFunctionBegin; 2679566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)pc,&comm)); 2689566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm,&size)); 269ea91fabdSFande Kong /* we do not need a hierarchical partitioning when 270ea91fabdSFande Kong * the total number of subdomains is consistent with 271ea91fabdSFande Kong * the number of MPI tasks. 272ea91fabdSFande Kong * For the following cases, we do not need to use HP 273ea91fabdSFande Kong * */ 274670417b2SFande Kong if (osm->N==PETSC_DETERMINE || osm->N>=size || osm->N==1) PetscFunctionReturn(0); 27563a3b9bcSJacob 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); 276ea91fabdSFande Kong nlocalsubdomains = size/osm->N; 277ea91fabdSFande Kong osm->n = 1; 2789566063dSJacob Faibussowitsch PetscCall(MatPartitioningCreate(comm,&part)); 2799566063dSJacob Faibussowitsch PetscCall(MatPartitioningSetAdjacency(part,pc->pmat)); 2809566063dSJacob Faibussowitsch PetscCall(MatPartitioningSetType(part,MATPARTITIONINGHIERARCH)); 2819566063dSJacob Faibussowitsch PetscCall(MatPartitioningHierarchicalSetNcoarseparts(part,osm->N)); 2829566063dSJacob Faibussowitsch PetscCall(MatPartitioningHierarchicalSetNfineparts(part,nlocalsubdomains)); 2839566063dSJacob Faibussowitsch PetscCall(MatPartitioningSetFromOptions(part)); 284ea91fabdSFande Kong /* get new processor owner number of each vertex */ 2859566063dSJacob Faibussowitsch PetscCall(MatPartitioningApply(part,&partitioning)); 2869566063dSJacob Faibussowitsch PetscCall(ISBuildTwoSided(partitioning,NULL,&fromrows)); 2879566063dSJacob Faibussowitsch PetscCall(ISPartitioningToNumbering(partitioning,&isn)); 2889566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isn)); 2899566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(fromrows,&fromrows_localsize)); 2909566063dSJacob Faibussowitsch PetscCall(MatPartitioningDestroy(&part)); 2919566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(pc->pmat,&outervec,NULL)); 2929566063dSJacob Faibussowitsch PetscCall(VecCreateMPI(comm,fromrows_localsize,PETSC_DETERMINE,&(osm->pcx))); 2939566063dSJacob Faibussowitsch PetscCall(VecDuplicate(osm->pcx,&(osm->pcy))); 2949566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(osm->pcx,NULL,outervec,fromrows,&(osm->pctoouter))); 2959566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(pc->pmat,fromrows,fromrows,MAT_INITIAL_MATRIX,&(osm->permutationP))); 2969566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)fromrows)); 297ea91fabdSFande Kong osm->permutationIS = fromrows; 298ea91fabdSFande Kong osm->pcmat = pc->pmat; 2999566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)osm->permutationP)); 300ea91fabdSFande Kong pc->pmat = osm->permutationP; 3019566063dSJacob Faibussowitsch PetscCall(VecDestroy(&outervec)); 3029566063dSJacob Faibussowitsch PetscCall(ISDestroy(&fromrows)); 3039566063dSJacob Faibussowitsch PetscCall(ISDestroy(&partitioning)); 304ea91fabdSFande Kong osm->n = PETSC_DETERMINE; 305ea91fabdSFande Kong PetscFunctionReturn(0); 306ea91fabdSFande Kong } 307ea91fabdSFande Kong 308f746d493SDmitry Karpeev static PetscErrorCode PCSetUp_GASM(PC pc) 309b1a0a8a3SJed Brown { 310f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 311ea91fabdSFande Kong PetscInt i,nInnerIndices,nTotalInnerIndices; 312c10bc1a1SDmitry Karpeev PetscMPIInt rank, size; 313b1a0a8a3SJed Brown MatReuse scall = MAT_REUSE_MATRIX; 314b1a0a8a3SJed Brown KSP ksp; 315b1a0a8a3SJed Brown PC subpc; 316b1a0a8a3SJed Brown const char *prefix,*pprefix; 317f746d493SDmitry Karpeev Vec x,y; 3186a4f0f73SDmitry Karpeev PetscInt oni; /* Number of indices in the i-th local outer subdomain. */ 3196a4f0f73SDmitry Karpeev const PetscInt *oidxi; /* Indices from the i-th subdomain local outer subdomain. */ 3206a4f0f73SDmitry Karpeev PetscInt on; /* Number of indices in the disjoint union of local outer subdomains. */ 3216a4f0f73SDmitry Karpeev PetscInt *oidx; /* Indices in the disjoint union of local outer subdomains. */ 3226a4f0f73SDmitry Karpeev IS gois; /* Disjoint union the global indices of outer subdomains. */ 3236a4f0f73SDmitry Karpeev IS goid; /* Identity IS of the size of the disjoint union of outer subdomains. */ 324f746d493SDmitry Karpeev PetscScalar *gxarray, *gyarray; 325930d09c1SFande Kong PetscInt gostart; /* Start of locally-owned indices in the vectors -- osm->gx,osm->gy -- over the disjoint union of outer subdomains. */ 3268f3b4b4dSDmitry Karpeev PetscInt num_subdomains = 0; 3270298fd71SBarry Smith DM *subdomain_dm = NULL; 3288f3b4b4dSDmitry Karpeev char **subdomain_names = NULL; 329f771a274SFande Kong PetscInt *numbering; 3308f3b4b4dSDmitry Karpeev 331b1a0a8a3SJed Brown PetscFunctionBegin; 3329566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size)); 3339566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank)); 334b1a0a8a3SJed Brown if (!pc->setupcalled) { 335ea91fabdSFande Kong /* use a hierarchical partitioning */ 336ea91fabdSFande Kong if (osm->hierarchicalpartitioning) { 3379566063dSJacob Faibussowitsch PetscCall(PCGASMSetHierarchicalPartitioning(pc)); 338ea91fabdSFande Kong } 3398f3b4b4dSDmitry Karpeev if (osm->n == PETSC_DETERMINE) { 3408f3b4b4dSDmitry Karpeev if (osm->N != PETSC_DETERMINE) { 3418f3b4b4dSDmitry Karpeev /* No local subdomains given, but the desired number of total subdomains is known, so construct them accordingly. */ 3429566063dSJacob Faibussowitsch PetscCall(PCGASMCreateSubdomains(pc->pmat,osm->N,&osm->n,&osm->iis)); 3438f3b4b4dSDmitry Karpeev } else if (osm->dm_subdomains && pc->dm) { 3448f3b4b4dSDmitry Karpeev /* try pc->dm next, if allowed */ 3458f3b4b4dSDmitry Karpeev PetscInt d; 346a35b7b57SDmitry Karpeev IS *inner_subdomain_is, *outer_subdomain_is; 3479566063dSJacob Faibussowitsch PetscCall(DMCreateDomainDecomposition(pc->dm, &num_subdomains, &subdomain_names, &inner_subdomain_is, &outer_subdomain_is, &subdomain_dm)); 348a35b7b57SDmitry Karpeev if (num_subdomains) { 3499566063dSJacob Faibussowitsch PetscCall(PCGASMSetSubdomains(pc, num_subdomains, inner_subdomain_is, outer_subdomain_is)); 35069ca1f37SDmitry Karpeev } 351a35b7b57SDmitry Karpeev for (d = 0; d < num_subdomains; ++d) { 3529566063dSJacob Faibussowitsch if (inner_subdomain_is) PetscCall(ISDestroy(&inner_subdomain_is[d])); 3539566063dSJacob Faibussowitsch if (outer_subdomain_is) PetscCall(ISDestroy(&outer_subdomain_is[d])); 354fdc48646SDmitry Karpeev } 3559566063dSJacob Faibussowitsch PetscCall(PetscFree(inner_subdomain_is)); 3569566063dSJacob Faibussowitsch PetscCall(PetscFree(outer_subdomain_is)); 3578f3b4b4dSDmitry Karpeev } else { 3588f3b4b4dSDmitry Karpeev /* still no subdomains; use one per processor */ 35944fe18e5SDmitry Karpeev osm->nmax = osm->n = 1; 3609566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size)); 361f746d493SDmitry Karpeev osm->N = size; 3629566063dSJacob Faibussowitsch PetscCall(PCGASMCreateLocalSubdomains(pc->pmat,osm->n,&osm->iis)); 363fdc48646SDmitry Karpeev } 36406b43e7eSDmitry Karpeev } 365a35b7b57SDmitry Karpeev if (!osm->iis) { 366a35b7b57SDmitry Karpeev /* 3678f3b4b4dSDmitry Karpeev osm->n was set in PCGASMSetSubdomains(), but the actual subdomains have not been supplied. 3688f3b4b4dSDmitry Karpeev We create the requisite number of local inner subdomains and then expand them into 3698f3b4b4dSDmitry Karpeev out subdomains, if necessary. 370a35b7b57SDmitry Karpeev */ 3719566063dSJacob Faibussowitsch PetscCall(PCGASMCreateLocalSubdomains(pc->pmat,osm->n,&osm->iis)); 372a35b7b57SDmitry Karpeev } 3738f3b4b4dSDmitry Karpeev if (!osm->ois) { 3748f3b4b4dSDmitry Karpeev /* 3758f3b4b4dSDmitry Karpeev Initially make outer subdomains the same as inner subdomains. If nonzero additional overlap 3768f3b4b4dSDmitry Karpeev has been requested, copy the inner subdomains over so they can be modified. 3778f3b4b4dSDmitry Karpeev */ 3789566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(osm->n,&osm->ois)); 3798f3b4b4dSDmitry Karpeev for (i=0; i<osm->n; ++i) { 380ea91fabdSFande Kong if (osm->overlap > 0 && osm->N>1) { /* With positive overlap, osm->iis[i] will be modified */ 3819566063dSJacob Faibussowitsch PetscCall(ISDuplicate(osm->iis[i],(osm->ois)+i)); 3829566063dSJacob Faibussowitsch PetscCall(ISCopy(osm->iis[i],osm->ois[i])); 3838f3b4b4dSDmitry Karpeev } else { 3849566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)((osm->iis)[i]))); 3858f3b4b4dSDmitry Karpeev osm->ois[i] = osm->iis[i]; 3868f3b4b4dSDmitry Karpeev } 3878f3b4b4dSDmitry Karpeev } 388ea91fabdSFande Kong if (osm->overlap>0 && osm->N>1) { 3898f3b4b4dSDmitry Karpeev /* Extend the "overlapping" regions by a number of steps */ 3909566063dSJacob Faibussowitsch PetscCall(MatIncreaseOverlapSplit(pc->pmat,osm->n,osm->ois,osm->overlap)); 3918f3b4b4dSDmitry Karpeev } 392b1a0a8a3SJed Brown } 3936a4f0f73SDmitry Karpeev 3948f3b4b4dSDmitry Karpeev /* Now the subdomains are defined. Determine their global and max local numbers, if necessary. */ 3958f3b4b4dSDmitry Karpeev if (osm->nmax == PETSC_DETERMINE) { 3968f3b4b4dSDmitry Karpeev PetscMPIInt inwork,outwork; 3978f3b4b4dSDmitry Karpeev /* determine global number of subdomains and the max number of local subdomains */ 3988f3b4b4dSDmitry Karpeev inwork = osm->n; 3991c2dc1cbSBarry Smith PetscCall(MPIU_Allreduce(&inwork,&outwork,1,MPI_INT,MPI_MAX,PetscObjectComm((PetscObject)pc))); 4008f3b4b4dSDmitry Karpeev osm->nmax = outwork; 4018f3b4b4dSDmitry Karpeev } 4028f3b4b4dSDmitry Karpeev if (osm->N == PETSC_DETERMINE) { 4038f3b4b4dSDmitry Karpeev /* Determine the number of globally-distinct subdomains and compute a global numbering for them. */ 4049566063dSJacob Faibussowitsch PetscCall(PetscObjectsListGetGlobalNumbering(PetscObjectComm((PetscObject)pc),osm->n,(PetscObject*)osm->ois,&osm->N,NULL)); 4058f3b4b4dSDmitry Karpeev } 4068f3b4b4dSDmitry Karpeev 407b1a0a8a3SJed Brown if (osm->sort_indices) { 408f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 4099566063dSJacob Faibussowitsch PetscCall(ISSort(osm->ois[i])); 4109566063dSJacob Faibussowitsch PetscCall(ISSort(osm->iis[i])); 411b1a0a8a3SJed Brown } 412b1a0a8a3SJed Brown } 4139566063dSJacob Faibussowitsch PetscCall(PCGetOptionsPrefix(pc,&prefix)); 4149566063dSJacob Faibussowitsch PetscCall(PCGASMPrintSubdomains(pc)); 4158f3b4b4dSDmitry Karpeev 4166a4f0f73SDmitry Karpeev /* 4176a4f0f73SDmitry Karpeev Merge the ISs, create merged vectors and restrictions. 4186a4f0f73SDmitry Karpeev */ 4196a4f0f73SDmitry Karpeev /* Merge outer subdomain ISs and construct a restriction onto the disjoint union of local outer subdomains. */ 4206a4f0f73SDmitry Karpeev on = 0; 421f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 4229566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(osm->ois[i],&oni)); 4236a4f0f73SDmitry Karpeev on += oni; 424f746d493SDmitry Karpeev } 4259566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(on, &oidx)); 4266a4f0f73SDmitry Karpeev on = 0; 427930d09c1SFande Kong /* Merge local indices together */ 428f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 4299566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(osm->ois[i],&oni)); 4309566063dSJacob Faibussowitsch PetscCall(ISGetIndices(osm->ois[i],&oidxi)); 4319566063dSJacob Faibussowitsch PetscCall(PetscArraycpy(oidx+on,oidxi,oni)); 4329566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(osm->ois[i],&oidxi)); 4336a4f0f73SDmitry Karpeev on += oni; 434f746d493SDmitry Karpeev } 4359566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(((PetscObject)(pc))->comm,on,oidx,PETSC_OWN_POINTER,&gois)); 436ea91fabdSFande Kong nTotalInnerIndices = 0; 437ea91fabdSFande Kong for (i=0; i<osm->n; i++) { 4389566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(osm->iis[i],&nInnerIndices)); 439ea91fabdSFande Kong nTotalInnerIndices += nInnerIndices; 440ea91fabdSFande Kong } 4419566063dSJacob Faibussowitsch PetscCall(VecCreateMPI(((PetscObject)(pc))->comm,nTotalInnerIndices,PETSC_DETERMINE,&x)); 4429566063dSJacob Faibussowitsch PetscCall(VecDuplicate(x,&y)); 443ea91fabdSFande Kong 4449566063dSJacob Faibussowitsch PetscCall(VecCreateMPI(PetscObjectComm((PetscObject)pc),on,PETSC_DECIDE,&osm->gx)); 4459566063dSJacob Faibussowitsch PetscCall(VecDuplicate(osm->gx,&osm->gy)); 4469566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(osm->gx, &gostart, NULL)); 4479566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PetscObjectComm((PetscObject)pc),on,gostart,1, &goid)); 448930d09c1SFande Kong /* gois might indices not on local */ 4499566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(x,gois,osm->gx,goid, &(osm->gorestriction))); 4509566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(osm->n,&numbering)); 4519566063dSJacob Faibussowitsch PetscCall(PetscObjectsListGetGlobalNumbering(PetscObjectComm((PetscObject)pc),osm->n,(PetscObject*)osm->ois,NULL,numbering)); 4529566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 4539566063dSJacob Faibussowitsch PetscCall(ISDestroy(&gois)); 4542fa5cd67SKarl Rupp 4556a4f0f73SDmitry Karpeev /* Merge inner subdomain ISs and construct a restriction onto the disjoint union of local inner subdomains. */ 4562fa5cd67SKarl Rupp { 4572fa5cd67SKarl Rupp PetscInt ini; /* Number of indices the i-th a local inner subdomain. */ 4588966356dSPierre Jolivet PetscInt in; /* Number of indices in the disjoint union of local inner subdomains. */ 4596a4f0f73SDmitry Karpeev PetscInt *iidx; /* Global indices in the merged local inner subdomain. */ 4606a4f0f73SDmitry Karpeev PetscInt *ioidx; /* Global indices of the disjoint union of inner subdomains within the disjoint union of outer subdomains. */ 4616a4f0f73SDmitry Karpeev IS giis; /* IS for the disjoint union of inner subdomains. */ 4626a4f0f73SDmitry Karpeev IS giois; /* IS for the disjoint union of inner subdomains within the disjoint union of outer subdomains. */ 463f771a274SFande Kong PetscScalar *array; 464f771a274SFande Kong const PetscInt *indices; 465f771a274SFande Kong PetscInt k; 4666a4f0f73SDmitry Karpeev on = 0; 467f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 4689566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(osm->ois[i],&oni)); 4696a4f0f73SDmitry Karpeev on += oni; 470f746d493SDmitry Karpeev } 4719566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(on, &iidx)); 4729566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(on, &ioidx)); 4739566063dSJacob Faibussowitsch PetscCall(VecGetArray(y,&array)); 474e12b4729SFande Kong /* set communicator id to determine where overlap is */ 475f771a274SFande Kong in = 0; 476f771a274SFande Kong for (i=0; i<osm->n; i++) { 4779566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(osm->iis[i],&ini)); 478f771a274SFande Kong for (k = 0; k < ini; ++k) { 479f771a274SFande Kong array[in+k] = numbering[i]; 480f771a274SFande Kong } 481f771a274SFande Kong in += ini; 482f771a274SFande Kong } 4839566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(y,&array)); 4849566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->gorestriction,y,osm->gy,INSERT_VALUES,SCATTER_FORWARD)); 4859566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->gorestriction,y,osm->gy,INSERT_VALUES,SCATTER_FORWARD)); 4869566063dSJacob Faibussowitsch PetscCall(VecGetOwnershipRange(osm->gy,&gostart, NULL)); 4879566063dSJacob Faibussowitsch PetscCall(VecGetArray(osm->gy,&array)); 488f771a274SFande Kong on = 0; 489f771a274SFande Kong in = 0; 490f771a274SFande Kong for (i=0; i<osm->n; i++) { 4919566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(osm->ois[i],&oni)); 4929566063dSJacob Faibussowitsch PetscCall(ISGetIndices(osm->ois[i],&indices)); 493f771a274SFande Kong for (k=0; k<oni; k++) { 494e12b4729SFande Kong /* skip overlapping indices to get inner domain */ 49543081b6cSDmitry Karpeev if (PetscRealPart(array[on+k]) != numbering[i]) continue; 496f771a274SFande Kong iidx[in] = indices[k]; 497f771a274SFande Kong ioidx[in++] = gostart+on+k; 498f771a274SFande Kong } 4999566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(osm->ois[i], &indices)); 500f771a274SFande Kong on += oni; 501f771a274SFande Kong } 5029566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(osm->gy,&array)); 5039566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)pc),in,iidx,PETSC_OWN_POINTER,&giis)); 5049566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PetscObjectComm((PetscObject)pc),in,ioidx,PETSC_OWN_POINTER,&giois)); 5059566063dSJacob Faibussowitsch PetscCall(VecScatterCreate(y,giis,osm->gy,giois,&osm->girestriction)); 5069566063dSJacob Faibussowitsch PetscCall(VecDestroy(&y)); 5079566063dSJacob Faibussowitsch PetscCall(ISDestroy(&giis)); 5089566063dSJacob Faibussowitsch PetscCall(ISDestroy(&giois)); 509b1a0a8a3SJed Brown } 5109566063dSJacob Faibussowitsch PetscCall(ISDestroy(&goid)); 5119566063dSJacob Faibussowitsch PetscCall(PetscFree(numbering)); 5122fa5cd67SKarl Rupp 513f746d493SDmitry Karpeev /* Create the subdomain work vectors. */ 5149566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(osm->n,&osm->x)); 5159566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(osm->n,&osm->y)); 5169566063dSJacob Faibussowitsch PetscCall(VecGetArray(osm->gx, &gxarray)); 5179566063dSJacob Faibussowitsch PetscCall(VecGetArray(osm->gy, &gyarray)); 5186a4f0f73SDmitry Karpeev for (i=0, on=0; i<osm->n; ++i, on += oni) { 5196a4f0f73SDmitry Karpeev PetscInt oNi; 5209566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(osm->ois[i],&oni)); 521930d09c1SFande Kong /* on a sub communicator */ 5229566063dSJacob Faibussowitsch PetscCall(ISGetSize(osm->ois[i],&oNi)); 5239566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(((PetscObject)(osm->ois[i]))->comm,1,oni,oNi,gxarray+on,&osm->x[i])); 5249566063dSJacob Faibussowitsch PetscCall(VecCreateMPIWithArray(((PetscObject)(osm->ois[i]))->comm,1,oni,oNi,gyarray+on,&osm->y[i])); 525b1a0a8a3SJed Brown } 5269566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(osm->gx, &gxarray)); 5279566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(osm->gy, &gyarray)); 5288f3b4b4dSDmitry Karpeev /* Create the subdomain solvers */ 5299566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(osm->n,&osm->ksp)); 5308f3b4b4dSDmitry Karpeev for (i=0; i<osm->n; i++) { 5318f3b4b4dSDmitry Karpeev char subprefix[PETSC_MAX_PATH_LEN+1]; 5329566063dSJacob Faibussowitsch PetscCall(KSPCreate(((PetscObject)(osm->ois[i]))->comm,&ksp)); 5339566063dSJacob Faibussowitsch PetscCall(KSPSetErrorIfNotConverged(ksp,pc->erroriffailure)); 5349566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)ksp)); 5359566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,1)); 5369566063dSJacob Faibussowitsch PetscCall(KSPSetType(ksp,KSPPREONLY)); 5379566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp,&subpc)); /* Why do we need this here? */ 5388f3b4b4dSDmitry Karpeev if (subdomain_dm) { 5399566063dSJacob Faibussowitsch PetscCall(KSPSetDM(ksp,subdomain_dm[i])); 5409566063dSJacob Faibussowitsch PetscCall(DMDestroy(subdomain_dm+i)); 5418f3b4b4dSDmitry Karpeev } 5429566063dSJacob Faibussowitsch PetscCall(PCGetOptionsPrefix(pc,&prefix)); 5439566063dSJacob Faibussowitsch PetscCall(KSPSetOptionsPrefix(ksp,prefix)); 5448f3b4b4dSDmitry Karpeev if (subdomain_names && subdomain_names[i]) { 5459566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(subprefix,PETSC_MAX_PATH_LEN,"sub_%s_",subdomain_names[i])); 5469566063dSJacob Faibussowitsch PetscCall(KSPAppendOptionsPrefix(ksp,subprefix)); 5479566063dSJacob Faibussowitsch PetscCall(PetscFree(subdomain_names[i])); 5488f3b4b4dSDmitry Karpeev } 5499566063dSJacob Faibussowitsch PetscCall(KSPAppendOptionsPrefix(ksp,"sub_")); 550b1a0a8a3SJed Brown osm->ksp[i] = ksp; 551b1a0a8a3SJed Brown } 5529566063dSJacob Faibussowitsch PetscCall(PetscFree(subdomain_dm)); 5539566063dSJacob Faibussowitsch PetscCall(PetscFree(subdomain_names)); 554b1a0a8a3SJed Brown scall = MAT_INITIAL_MATRIX; 5558f3b4b4dSDmitry Karpeev } else { /* if (pc->setupcalled) */ 556b1a0a8a3SJed Brown /* 5578f3b4b4dSDmitry Karpeev Destroy the submatrices from the previous iteration 558b1a0a8a3SJed Brown */ 559b1a0a8a3SJed Brown if (pc->flag == DIFFERENT_NONZERO_PATTERN) { 5609566063dSJacob Faibussowitsch PetscCall(MatDestroyMatrices(osm->n,&osm->pmat)); 561b1a0a8a3SJed Brown scall = MAT_INITIAL_MATRIX; 562b1a0a8a3SJed Brown } 563ea91fabdSFande Kong if (osm->permutationIS) { 5649566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrix(pc->pmat,osm->permutationIS,osm->permutationIS,scall,&osm->permutationP)); 5659566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)osm->permutationP)); 566ea91fabdSFande Kong osm->pcmat = pc->pmat; 567ea91fabdSFande Kong pc->pmat = osm->permutationP; 568b1a0a8a3SJed Brown } 569ea91fabdSFande Kong } 570ea91fabdSFande Kong 571b1a0a8a3SJed Brown /* 5722da392ccSBarry Smith Extract the submatrices. 573b1a0a8a3SJed Brown */ 57481a5c4aaSDmitry Karpeev if (size > 1) { 5759566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatricesMPI(pc->pmat,osm->n,osm->ois,osm->ois,scall,&osm->pmat)); 5762fa5cd67SKarl Rupp } else { 5779566063dSJacob Faibussowitsch PetscCall(MatCreateSubMatrices(pc->pmat,osm->n,osm->ois,osm->ois,scall,&osm->pmat)); 57881a5c4aaSDmitry Karpeev } 579b1a0a8a3SJed Brown if (scall == MAT_INITIAL_MATRIX) { 5809566063dSJacob Faibussowitsch PetscCall(PetscObjectGetOptionsPrefix((PetscObject)pc->pmat,&pprefix)); 581f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 5829566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)pc,(PetscObject)osm->pmat[i])); 5839566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)osm->pmat[i],pprefix)); 584b1a0a8a3SJed Brown } 585b1a0a8a3SJed Brown } 586b1a0a8a3SJed Brown 587b1a0a8a3SJed Brown /* Return control to the user so that the submatrices can be modified (e.g., to apply 588b1a0a8a3SJed Brown different boundary conditions for the submatrices than for the global problem) */ 5899566063dSJacob Faibussowitsch PetscCall(PCModifySubMatrices(pc,osm->n,osm->ois,osm->ois,osm->pmat,pc->modifysubmatricesP)); 590b1a0a8a3SJed Brown 591b1a0a8a3SJed Brown /* 5926a4f0f73SDmitry Karpeev Loop over submatrices putting them into local ksps 593b1a0a8a3SJed Brown */ 594f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 5959566063dSJacob Faibussowitsch PetscCall(KSPSetOperators(osm->ksp[i],osm->pmat[i],osm->pmat[i])); 5969566063dSJacob Faibussowitsch PetscCall(KSPGetOptionsPrefix(osm->ksp[i],&prefix)); 5979566063dSJacob Faibussowitsch PetscCall(MatSetOptionsPrefix(osm->pmat[i],prefix)); 598b1a0a8a3SJed Brown if (!pc->setupcalled) { 5999566063dSJacob Faibussowitsch PetscCall(KSPSetFromOptions(osm->ksp[i])); 600b1a0a8a3SJed Brown } 601b1a0a8a3SJed Brown } 602ea91fabdSFande Kong if (osm->pcmat) { 6039566063dSJacob Faibussowitsch PetscCall(MatDestroy(&pc->pmat)); 604ea91fabdSFande Kong pc->pmat = osm->pcmat; 6050a545947SLisandro Dalcin osm->pcmat = NULL; 606ea91fabdSFande Kong } 607b1a0a8a3SJed Brown PetscFunctionReturn(0); 608b1a0a8a3SJed Brown } 609b1a0a8a3SJed Brown 610f746d493SDmitry Karpeev static PetscErrorCode PCSetUpOnBlocks_GASM(PC pc) 611b1a0a8a3SJed Brown { 612f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 613b1a0a8a3SJed Brown PetscInt i; 614b1a0a8a3SJed Brown 615b1a0a8a3SJed Brown PetscFunctionBegin; 616f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 6179566063dSJacob Faibussowitsch PetscCall(KSPSetUp(osm->ksp[i])); 618b1a0a8a3SJed Brown } 619b1a0a8a3SJed Brown PetscFunctionReturn(0); 620b1a0a8a3SJed Brown } 621b1a0a8a3SJed Brown 622ea91fabdSFande Kong static PetscErrorCode PCApply_GASM(PC pc,Vec xin,Vec yout) 623b1a0a8a3SJed Brown { 624f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 625f746d493SDmitry Karpeev PetscInt i; 626ea91fabdSFande Kong Vec x,y; 627b1a0a8a3SJed Brown ScatterMode forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE; 628b1a0a8a3SJed Brown 629b1a0a8a3SJed Brown PetscFunctionBegin; 630ea91fabdSFande Kong if (osm->pctoouter) { 6319566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->pctoouter,xin,osm->pcx,INSERT_VALUES,SCATTER_REVERSE)); 6329566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->pctoouter,xin,osm->pcx,INSERT_VALUES,SCATTER_REVERSE)); 633ea91fabdSFande Kong x = osm->pcx; 634ea91fabdSFande Kong y = osm->pcy; 635ea91fabdSFande Kong } else { 636ea91fabdSFande Kong x = xin; 637ea91fabdSFande Kong y = yout; 638ea91fabdSFande Kong } 639b1a0a8a3SJed Brown /* 64048e38a8aSPierre Jolivet support for limiting the restriction or interpolation only to the inner 641b1a0a8a3SJed Brown subdomain values (leaving the other values 0). 642b1a0a8a3SJed Brown */ 643f746d493SDmitry Karpeev if (!(osm->type & PC_GASM_RESTRICT)) { 644b1a0a8a3SJed Brown /* have to zero the work RHS since scatter may leave some slots empty */ 6459566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(osm->gx)); 6469566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->girestriction,x,osm->gx,INSERT_VALUES,forward)); 6472fa5cd67SKarl Rupp } else { 6489566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward)); 649b1a0a8a3SJed Brown } 6509566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(osm->gy)); 6516a4f0f73SDmitry Karpeev if (!(osm->type & PC_GASM_RESTRICT)) { 6529566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->girestriction,x,osm->gx,INSERT_VALUES,forward)); 6532fa5cd67SKarl Rupp } else { 6549566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward)); 6556a4f0f73SDmitry Karpeev } 65686b96d36SDmitry Karpeev /* do the subdomain solves */ 65786b96d36SDmitry Karpeev for (i=0; i<osm->n; ++i) { 6589566063dSJacob Faibussowitsch PetscCall(KSPSolve(osm->ksp[i],osm->x[i],osm->y[i])); 6599566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(osm->ksp[i],pc,osm->y[i])); 660b1a0a8a3SJed Brown } 66148e38a8aSPierre Jolivet /* do we need to zero y? */ 6629566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(y)); 6636a4f0f73SDmitry Karpeev if (!(osm->type & PC_GASM_INTERPOLATE)) { 6649566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->girestriction,osm->gy,y,ADD_VALUES,reverse)); 6659566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->girestriction,osm->gy,y,ADD_VALUES,reverse)); 6662fa5cd67SKarl Rupp } else { 6679566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse)); 6689566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse)); 6696a4f0f73SDmitry Karpeev } 670ea91fabdSFande Kong if (osm->pctoouter) { 6719566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->pctoouter,y,yout,INSERT_VALUES,SCATTER_FORWARD)); 6729566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->pctoouter,y,yout,INSERT_VALUES,SCATTER_FORWARD)); 673ea91fabdSFande Kong } 674ea91fabdSFande Kong PetscFunctionReturn(0); 675b1a0a8a3SJed Brown } 676b1a0a8a3SJed Brown 67748e38a8aSPierre Jolivet static PetscErrorCode PCMatApply_GASM(PC pc,Mat Xin,Mat Yout) 67848e38a8aSPierre Jolivet { 67948e38a8aSPierre Jolivet PC_GASM *osm = (PC_GASM*)pc->data; 68048e38a8aSPierre Jolivet Mat X,Y,O=NULL,Z,W; 68148e38a8aSPierre Jolivet Vec x,y; 68248e38a8aSPierre Jolivet PetscInt i,m,M,N; 68348e38a8aSPierre Jolivet ScatterMode forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE; 68448e38a8aSPierre Jolivet 68548e38a8aSPierre Jolivet PetscFunctionBegin; 68608401ef6SPierre Jolivet PetscCheck(osm->n == 1,PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Not yet implemented"); 6879566063dSJacob Faibussowitsch PetscCall(MatGetSize(Xin,NULL,&N)); 68848e38a8aSPierre Jolivet if (osm->pctoouter) { 6899566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(osm->pcx,&m)); 6909566063dSJacob Faibussowitsch PetscCall(VecGetSize(osm->pcx,&M)); 6919566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PetscObjectComm((PetscObject)osm->ois[0]),m,PETSC_DECIDE,M,N,NULL,&O)); 69248e38a8aSPierre Jolivet for (i = 0; i < N; ++i) { 6939566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumnVecRead(Xin,i,&x)); 6949566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumnVecWrite(O,i,&y)); 6959566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->pctoouter,x,y,INSERT_VALUES,SCATTER_REVERSE)); 6969566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->pctoouter,x,y,INSERT_VALUES,SCATTER_REVERSE)); 6979566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumnVecWrite(O,i,&y)); 6989566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumnVecRead(Xin,i,&x)); 69948e38a8aSPierre Jolivet } 70048e38a8aSPierre Jolivet X = Y = O; 70148e38a8aSPierre Jolivet } else { 70248e38a8aSPierre Jolivet X = Xin; 70348e38a8aSPierre Jolivet Y = Yout; 70448e38a8aSPierre Jolivet } 70548e38a8aSPierre Jolivet /* 70648e38a8aSPierre Jolivet support for limiting the restriction or interpolation only to the inner 70748e38a8aSPierre Jolivet subdomain values (leaving the other values 0). 70848e38a8aSPierre Jolivet */ 7099566063dSJacob Faibussowitsch PetscCall(VecGetLocalSize(osm->x[0],&m)); 7109566063dSJacob Faibussowitsch PetscCall(VecGetSize(osm->x[0],&M)); 7119566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PetscObjectComm((PetscObject)osm->ois[0]),m,PETSC_DECIDE,M,N,NULL,&Z)); 71248e38a8aSPierre Jolivet for (i = 0; i < N; ++i) { 7139566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumnVecRead(X,i,&x)); 7149566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumnVecWrite(Z,i,&y)); 71548e38a8aSPierre Jolivet if (!(osm->type & PC_GASM_RESTRICT)) { 71648e38a8aSPierre Jolivet /* have to zero the work RHS since scatter may leave some slots empty */ 7179566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(y)); 7189566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->girestriction,x,y,INSERT_VALUES,forward)); 7199566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->girestriction,x,y,INSERT_VALUES,forward)); 72048e38a8aSPierre Jolivet } else { 7219566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->gorestriction,x,y,INSERT_VALUES,forward)); 7229566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->gorestriction,x,y,INSERT_VALUES,forward)); 72348e38a8aSPierre Jolivet } 7249566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumnVecWrite(Z,i,&y)); 7259566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumnVecRead(X,i,&x)); 72648e38a8aSPierre Jolivet } 7279566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PetscObjectComm((PetscObject)osm->ois[0]),m,PETSC_DECIDE,M,N,NULL,&W)); 7289566063dSJacob Faibussowitsch PetscCall(MatSetOption(Z,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE)); 7299566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(Z,MAT_FINAL_ASSEMBLY)); 7309566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(Z,MAT_FINAL_ASSEMBLY)); 73148e38a8aSPierre Jolivet /* do the subdomain solve */ 7329566063dSJacob Faibussowitsch PetscCall(KSPMatSolve(osm->ksp[0],Z,W)); 7339566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(osm->ksp[0],pc,NULL)); 7349566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Z)); 73548e38a8aSPierre Jolivet /* do we need to zero y? */ 7369566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(Y)); 73748e38a8aSPierre Jolivet for (i = 0; i < N; ++i) { 7389566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumnVecWrite(Y,i,&y)); 7399566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumnVecRead(W,i,&x)); 74048e38a8aSPierre Jolivet if (!(osm->type & PC_GASM_INTERPOLATE)) { 7419566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->girestriction,x,y,ADD_VALUES,reverse)); 7429566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->girestriction,x,y,ADD_VALUES,reverse)); 74348e38a8aSPierre Jolivet } else { 7449566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->gorestriction,x,y,ADD_VALUES,reverse)); 7459566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->gorestriction,x,y,ADD_VALUES,reverse)); 74648e38a8aSPierre Jolivet } 7479566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumnVecRead(W,i,&x)); 74848e38a8aSPierre Jolivet if (osm->pctoouter) { 7499566063dSJacob Faibussowitsch PetscCall(MatDenseGetColumnVecWrite(Yout,i,&x)); 7509566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->pctoouter,y,x,INSERT_VALUES,SCATTER_FORWARD)); 7519566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->pctoouter,y,x,INSERT_VALUES,SCATTER_FORWARD)); 7529566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumnVecRead(Yout,i,&x)); 75348e38a8aSPierre Jolivet } 7549566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreColumnVecWrite(Y,i,&y)); 75548e38a8aSPierre Jolivet } 7569566063dSJacob Faibussowitsch PetscCall(MatDestroy(&W)); 7579566063dSJacob Faibussowitsch PetscCall(MatDestroy(&O)); 75848e38a8aSPierre Jolivet PetscFunctionReturn(0); 75948e38a8aSPierre Jolivet } 76048e38a8aSPierre Jolivet 761ea91fabdSFande Kong static PetscErrorCode PCApplyTranspose_GASM(PC pc,Vec xin,Vec yout) 762b1a0a8a3SJed Brown { 763f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 764f746d493SDmitry Karpeev PetscInt i; 765ea91fabdSFande Kong Vec x,y; 766b1a0a8a3SJed Brown ScatterMode forward = SCATTER_FORWARD,reverse = SCATTER_REVERSE; 767b1a0a8a3SJed Brown 768b1a0a8a3SJed Brown PetscFunctionBegin; 769ea91fabdSFande Kong if (osm->pctoouter) { 7709566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->pctoouter,xin,osm->pcx,INSERT_VALUES,SCATTER_REVERSE)); 7719566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->pctoouter,xin,osm->pcx,INSERT_VALUES,SCATTER_REVERSE)); 772ea91fabdSFande Kong x = osm->pcx; 773ea91fabdSFande Kong y = osm->pcy; 774ea91fabdSFande Kong }else{ 775ea91fabdSFande Kong x = xin; 776ea91fabdSFande Kong y = yout; 777ea91fabdSFande Kong } 778b1a0a8a3SJed Brown /* 779b1a0a8a3SJed Brown Support for limiting the restriction or interpolation to only local 780b1a0a8a3SJed Brown subdomain values (leaving the other values 0). 781b1a0a8a3SJed Brown 782f746d493SDmitry Karpeev Note: these are reversed from the PCApply_GASM() because we are applying the 783b1a0a8a3SJed Brown transpose of the three terms 784b1a0a8a3SJed Brown */ 785f746d493SDmitry Karpeev if (!(osm->type & PC_GASM_INTERPOLATE)) { 786b1a0a8a3SJed Brown /* have to zero the work RHS since scatter may leave some slots empty */ 7879566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(osm->gx)); 7889566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->girestriction,x,osm->gx,INSERT_VALUES,forward)); 7892fa5cd67SKarl Rupp } else { 7909566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward)); 791b1a0a8a3SJed Brown } 7929566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(osm->gy)); 7936a4f0f73SDmitry Karpeev if (!(osm->type & PC_GASM_INTERPOLATE)) { 7949566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->girestriction,x,osm->gx,INSERT_VALUES,forward)); 7952fa5cd67SKarl Rupp } else { 7969566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward)); 7976a4f0f73SDmitry Karpeev } 798b1a0a8a3SJed Brown /* do the local solves */ 799ab3e923fSDmitry 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. */ 8009566063dSJacob Faibussowitsch PetscCall(KSPSolveTranspose(osm->ksp[i],osm->x[i],osm->y[i])); 8019566063dSJacob Faibussowitsch PetscCall(KSPCheckSolve(osm->ksp[i],pc,osm->y[i])); 802b1a0a8a3SJed Brown } 8039566063dSJacob Faibussowitsch PetscCall(VecZeroEntries(y)); 8046a4f0f73SDmitry Karpeev if (!(osm->type & PC_GASM_RESTRICT)) { 8059566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->girestriction,osm->gy,y,ADD_VALUES,reverse)); 8069566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->girestriction,osm->gy,y,ADD_VALUES,reverse)); 8072fa5cd67SKarl Rupp } else { 8089566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse)); 8099566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse)); 8106a4f0f73SDmitry Karpeev } 811ea91fabdSFande Kong if (osm->pctoouter) { 8129566063dSJacob Faibussowitsch PetscCall(VecScatterBegin(osm->pctoouter,y,yout,INSERT_VALUES,SCATTER_FORWARD)); 8139566063dSJacob Faibussowitsch PetscCall(VecScatterEnd(osm->pctoouter,y,yout,INSERT_VALUES,SCATTER_FORWARD)); 814ea91fabdSFande Kong } 815b1a0a8a3SJed Brown PetscFunctionReturn(0); 816b1a0a8a3SJed Brown } 817b1a0a8a3SJed Brown 818a06653b4SBarry Smith static PetscErrorCode PCReset_GASM(PC pc) 819b1a0a8a3SJed Brown { 820f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 821b1a0a8a3SJed Brown PetscInt i; 822b1a0a8a3SJed Brown 823b1a0a8a3SJed Brown PetscFunctionBegin; 824b1a0a8a3SJed Brown if (osm->ksp) { 825f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 8269566063dSJacob Faibussowitsch PetscCall(KSPReset(osm->ksp[i])); 827b1a0a8a3SJed Brown } 828b1a0a8a3SJed Brown } 829b1a0a8a3SJed Brown if (osm->pmat) { 830f746d493SDmitry Karpeev if (osm->n > 0) { 831df750dc8SHong Zhang PetscMPIInt size; 8329566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size)); 833df750dc8SHong Zhang if (size > 1) { 8347dae84e0SHong Zhang /* osm->pmat is created by MatCreateSubMatricesMPI(), cannot use MatDestroySubMatrices() */ 8359566063dSJacob Faibussowitsch PetscCall(MatDestroyMatrices(osm->n,&osm->pmat)); 836df750dc8SHong Zhang } else { 8379566063dSJacob Faibussowitsch PetscCall(MatDestroySubMatrices(osm->n,&osm->pmat)); 838df750dc8SHong Zhang } 839b1a0a8a3SJed Brown } 840b1a0a8a3SJed Brown } 841d34fcf5fSBarry Smith if (osm->x) { 842f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 8439566063dSJacob Faibussowitsch PetscCall(VecDestroy(&osm->x[i])); 8449566063dSJacob Faibussowitsch PetscCall(VecDestroy(&osm->y[i])); 845b1a0a8a3SJed Brown } 846d34fcf5fSBarry Smith } 8479566063dSJacob Faibussowitsch PetscCall(VecDestroy(&osm->gx)); 8489566063dSJacob Faibussowitsch PetscCall(VecDestroy(&osm->gy)); 849ab3e923fSDmitry Karpeev 8509566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&osm->gorestriction)); 8519566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&osm->girestriction)); 8528f3b4b4dSDmitry Karpeev if (!osm->user_subdomains) { 8539566063dSJacob Faibussowitsch PetscCall(PCGASMDestroySubdomains(osm->n,&osm->ois,&osm->iis)); 8548f3b4b4dSDmitry Karpeev osm->N = PETSC_DETERMINE; 8558f3b4b4dSDmitry Karpeev osm->nmax = PETSC_DETERMINE; 8568f3b4b4dSDmitry Karpeev } 857ea91fabdSFande Kong if (osm->pctoouter) { 8589566063dSJacob Faibussowitsch PetscCall(VecScatterDestroy(&(osm->pctoouter))); 859ea91fabdSFande Kong } 860ea91fabdSFande Kong if (osm->permutationIS) { 8619566063dSJacob Faibussowitsch PetscCall(ISDestroy(&(osm->permutationIS))); 862ea91fabdSFande Kong } 863ea91fabdSFande Kong if (osm->pcx) { 8649566063dSJacob Faibussowitsch PetscCall(VecDestroy(&(osm->pcx))); 865ea91fabdSFande Kong } 866ea91fabdSFande Kong if (osm->pcy) { 8679566063dSJacob Faibussowitsch PetscCall(VecDestroy(&(osm->pcy))); 868ea91fabdSFande Kong } 869ea91fabdSFande Kong if (osm->permutationP) { 8709566063dSJacob Faibussowitsch PetscCall(MatDestroy(&(osm->permutationP))); 871ea91fabdSFande Kong } 872ea91fabdSFande Kong if (osm->pcmat) { 8739566063dSJacob Faibussowitsch PetscCall(MatDestroy(&osm->pcmat)); 874ea91fabdSFande Kong } 875a06653b4SBarry Smith PetscFunctionReturn(0); 876a06653b4SBarry Smith } 877a06653b4SBarry Smith 878a06653b4SBarry Smith static PetscErrorCode PCDestroy_GASM(PC pc) 879a06653b4SBarry Smith { 880a06653b4SBarry Smith PC_GASM *osm = (PC_GASM*)pc->data; 881a06653b4SBarry Smith PetscInt i; 882a06653b4SBarry Smith 883a06653b4SBarry Smith PetscFunctionBegin; 8849566063dSJacob Faibussowitsch PetscCall(PCReset_GASM(pc)); 885135757f6SDmitry Karpeev /* PCReset will not destroy subdomains, if user_subdomains is true. */ 8869566063dSJacob Faibussowitsch PetscCall(PCGASMDestroySubdomains(osm->n,&osm->ois,&osm->iis)); 887a06653b4SBarry Smith if (osm->ksp) { 888a06653b4SBarry Smith for (i=0; i<osm->n; i++) { 8899566063dSJacob Faibussowitsch PetscCall(KSPDestroy(&osm->ksp[i])); 890a06653b4SBarry Smith } 8919566063dSJacob Faibussowitsch PetscCall(PetscFree(osm->ksp)); 892a06653b4SBarry Smith } 8939566063dSJacob Faibussowitsch PetscCall(PetscFree(osm->x)); 8949566063dSJacob Faibussowitsch PetscCall(PetscFree(osm->y)); 895*2e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetSubdomains_C",NULL)); 896*2e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetOverlap_C",NULL)); 897*2e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetType_C",NULL)); 898*2e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetSortIndices_C",NULL)); 899*2e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGASMGetSubKSP_C",NULL)); 9009566063dSJacob Faibussowitsch PetscCall(PetscFree(pc->data)); 901b1a0a8a3SJed Brown PetscFunctionReturn(0); 902b1a0a8a3SJed Brown } 903b1a0a8a3SJed Brown 9044416b707SBarry Smith static PetscErrorCode PCSetFromOptions_GASM(PetscOptionItems *PetscOptionsObject,PC pc) 905a6dfd86eSKarl Rupp { 906f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 907b1a0a8a3SJed Brown PetscInt blocks,ovl; 90857501b6eSBarry Smith PetscBool flg; 909f746d493SDmitry Karpeev PCGASMType gasmtype; 910b1a0a8a3SJed Brown 911b1a0a8a3SJed Brown PetscFunctionBegin; 912d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject,"Generalized additive Schwarz options"); 9139566063dSJacob 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)); 9149566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_gasm_total_subdomains","Total number of subdomains across communicator","PCGASMSetTotalSubdomains",osm->N,&blocks,&flg)); 91565db9045SDmitry Karpeev if (flg) { 9169566063dSJacob Faibussowitsch PetscCall(PCGASMSetTotalSubdomains(pc,blocks)); 91765db9045SDmitry Karpeev } 9189566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-pc_gasm_overlap","Number of overlapping degrees of freedom","PCGASMSetOverlap",osm->overlap,&ovl,&flg)); 91965db9045SDmitry Karpeev if (flg) { 9209566063dSJacob Faibussowitsch PetscCall(PCGASMSetOverlap(pc,ovl)); 921d709ea83SDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 92265db9045SDmitry Karpeev } 923b1a0a8a3SJed Brown flg = PETSC_FALSE; 9249566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-pc_gasm_type","Type of restriction/extension","PCGASMSetType",PCGASMTypes,(PetscEnum)osm->type,(PetscEnum*)&gasmtype,&flg)); 9259566063dSJacob Faibussowitsch if (flg) PetscCall(PCGASMSetType(pc,gasmtype)); 9269566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-pc_gasm_use_hierachical_partitioning","use hierarchical partitioning",NULL,osm->hierarchicalpartitioning,&osm->hierarchicalpartitioning,&flg)); 927d0609cedSBarry Smith PetscOptionsHeadEnd(); 928b1a0a8a3SJed Brown PetscFunctionReturn(0); 929b1a0a8a3SJed Brown } 930b1a0a8a3SJed Brown 931b1a0a8a3SJed Brown /*------------------------------------------------------------------------------------*/ 932b1a0a8a3SJed Brown 9338f3b4b4dSDmitry Karpeev /*@ 9348f3b4b4dSDmitry Karpeev PCGASMSetTotalSubdomains - sets the total number of subdomains to use across the 9358f3b4b4dSDmitry Karpeev communicator. 9368f3b4b4dSDmitry Karpeev Logically collective on pc 9378f3b4b4dSDmitry Karpeev 9388f3b4b4dSDmitry Karpeev Input Parameters: 9398f3b4b4dSDmitry Karpeev + pc - the preconditioner 9408f3b4b4dSDmitry Karpeev - N - total number of subdomains 9418f3b4b4dSDmitry Karpeev 9428f3b4b4dSDmitry Karpeev Level: beginner 9438f3b4b4dSDmitry Karpeev 944db781477SPatrick Sanan .seealso: `PCGASMSetSubdomains()`, `PCGASMSetOverlap()` 945db781477SPatrick Sanan `PCGASMCreateSubdomains2D()` 9468f3b4b4dSDmitry Karpeev @*/ 9478f3b4b4dSDmitry Karpeev PetscErrorCode PCGASMSetTotalSubdomains(PC pc,PetscInt N) 9488f3b4b4dSDmitry Karpeev { 9498f3b4b4dSDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 9508f3b4b4dSDmitry Karpeev PetscMPIInt size,rank; 9518f3b4b4dSDmitry Karpeev 9528f3b4b4dSDmitry Karpeev PetscFunctionBegin; 95363a3b9bcSJacob 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); 95428b400f6SJacob Faibussowitsch PetscCheck(!pc->setupcalled,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetTotalSubdomains() should be called before calling PCSetUp()."); 9558f3b4b4dSDmitry Karpeev 9569566063dSJacob Faibussowitsch PetscCall(PCGASMDestroySubdomains(osm->n,&osm->iis,&osm->ois)); 9578f3b4b4dSDmitry Karpeev osm->ois = osm->iis = NULL; 9588f3b4b4dSDmitry Karpeev 9599566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size)); 9609566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank)); 9618f3b4b4dSDmitry Karpeev osm->N = N; 9628f3b4b4dSDmitry Karpeev osm->n = PETSC_DETERMINE; 9638f3b4b4dSDmitry Karpeev osm->nmax = PETSC_DETERMINE; 9648f3b4b4dSDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 9658f3b4b4dSDmitry Karpeev PetscFunctionReturn(0); 9668f3b4b4dSDmitry Karpeev } 9678f3b4b4dSDmitry Karpeev 968b541e6a4SDmitry Karpeev static PetscErrorCode PCGASMSetSubdomains_GASM(PC pc,PetscInt n,IS iis[],IS ois[]) 969b1a0a8a3SJed Brown { 970f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 971b1a0a8a3SJed Brown PetscInt i; 972b1a0a8a3SJed Brown 973b1a0a8a3SJed Brown PetscFunctionBegin; 97463a3b9bcSJacob Faibussowitsch PetscCheck(n >= 1,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Each process must have 1 or more subdomains, got n = %" PetscInt_FMT,n); 97528b400f6SJacob Faibussowitsch PetscCheck(!pc->setupcalled,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetSubdomains() should be called before calling PCSetUp()."); 976b1a0a8a3SJed Brown 9779566063dSJacob Faibussowitsch PetscCall(PCGASMDestroySubdomains(osm->n,&osm->iis,&osm->ois)); 9788f3b4b4dSDmitry Karpeev osm->iis = osm->ois = NULL; 9798f3b4b4dSDmitry Karpeev osm->n = n; 9808f3b4b4dSDmitry Karpeev osm->N = PETSC_DETERMINE; 9818f3b4b4dSDmitry Karpeev osm->nmax = PETSC_DETERMINE; 982a35b7b57SDmitry Karpeev if (ois) { 9839566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n,&osm->ois)); 9848f3b4b4dSDmitry Karpeev for (i=0; i<n; i++) { 9859566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)ois[i])); 9868f3b4b4dSDmitry Karpeev osm->ois[i] = ois[i]; 9878f3b4b4dSDmitry Karpeev } 9888f3b4b4dSDmitry Karpeev /* 9898f3b4b4dSDmitry Karpeev Since the user set the outer subdomains, even if nontrivial overlap was requested via PCGASMSetOverlap(), 9908f3b4b4dSDmitry Karpeev it will be ignored. To avoid confusion later on (e.g., when viewing the PC), the overlap size is set to -1. 9918f3b4b4dSDmitry Karpeev */ 992b1a0a8a3SJed Brown osm->overlap = -1; 993670417b2SFande Kong /* inner subdomains must be provided */ 99428b400f6SJacob Faibussowitsch PetscCheck(iis,PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"inner indices have to be provided "); 995670417b2SFande Kong }/* end if */ 996a35b7b57SDmitry Karpeev if (iis) { 9979566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(n,&osm->iis)); 9988f3b4b4dSDmitry Karpeev for (i=0; i<n; i++) { 9999566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)iis[i])); 10008f3b4b4dSDmitry Karpeev osm->iis[i] = iis[i]; 10018f3b4b4dSDmitry Karpeev } 1002a35b7b57SDmitry Karpeev if (!ois) { 1003390e1bf2SBarry Smith osm->ois = NULL; 1004670417b2SFande Kong /* if user does not provide outer indices, we will create the corresponding outer indices using osm->overlap =1 in PCSetUp_GASM */ 1005670417b2SFande Kong } 1006670417b2SFande Kong } 100776bd3646SJed Brown if (PetscDefined(USE_DEBUG)) { 1008670417b2SFande Kong PetscInt j,rstart,rend,*covered,lsize; 1009670417b2SFande Kong const PetscInt *indices; 1010670417b2SFande Kong /* check if the inner indices cover and only cover the local portion of the preconditioning matrix */ 10119566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(pc->pmat,&rstart,&rend)); 10129566063dSJacob Faibussowitsch PetscCall(PetscCalloc1(rend-rstart,&covered)); 1013670417b2SFande Kong /* check if the current processor owns indices from others */ 1014a35b7b57SDmitry Karpeev for (i=0; i<n; i++) { 10159566063dSJacob Faibussowitsch PetscCall(ISGetIndices(osm->iis[i],&indices)); 10169566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(osm->iis[i],&lsize)); 1017670417b2SFande Kong for (j=0; j<lsize; j++) { 10182472a847SBarry 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]); 10192472a847SBarry 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]); 102063a3b9bcSJacob Faibussowitsch covered[indices[j]-rstart] = 1; 1021a35b7b57SDmitry Karpeev } 10229566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(osm->iis[i],&indices)); 10238f3b4b4dSDmitry Karpeev } 1024670417b2SFande Kong /* check if we miss any indices */ 1025670417b2SFande Kong for (i=rstart; i<rend; i++) { 102663a3b9bcSJacob Faibussowitsch PetscCheck(covered[i-rstart],PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"local entity %" PetscInt_FMT " was not covered by inner subdomains",i); 1027a35b7b57SDmitry Karpeev } 10289566063dSJacob Faibussowitsch PetscCall(PetscFree(covered)); 1029a35b7b57SDmitry Karpeev } 10308f3b4b4dSDmitry Karpeev if (iis) osm->user_subdomains = PETSC_TRUE; 1031b1a0a8a3SJed Brown PetscFunctionReturn(0); 1032b1a0a8a3SJed Brown } 1033b1a0a8a3SJed Brown 1034f7a08781SBarry Smith static PetscErrorCode PCGASMSetOverlap_GASM(PC pc,PetscInt ovl) 1035b1a0a8a3SJed Brown { 1036f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 1037b1a0a8a3SJed Brown 1038b1a0a8a3SJed Brown PetscFunctionBegin; 103908401ef6SPierre Jolivet PetscCheck(ovl >= 0,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap value requested"); 10402472a847SBarry Smith PetscCheck(!pc->setupcalled || ovl == osm->overlap,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetOverlap() should be called before PCSetUp()."); 10412fa5cd67SKarl Rupp if (!pc->setupcalled) osm->overlap = ovl; 1042b1a0a8a3SJed Brown PetscFunctionReturn(0); 1043b1a0a8a3SJed Brown } 1044b1a0a8a3SJed Brown 1045f7a08781SBarry Smith static PetscErrorCode PCGASMSetType_GASM(PC pc,PCGASMType type) 1046b1a0a8a3SJed Brown { 1047f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 1048b1a0a8a3SJed Brown 1049b1a0a8a3SJed Brown PetscFunctionBegin; 1050b1a0a8a3SJed Brown osm->type = type; 1051b1a0a8a3SJed Brown osm->type_set = PETSC_TRUE; 1052b1a0a8a3SJed Brown PetscFunctionReturn(0); 1053b1a0a8a3SJed Brown } 1054b1a0a8a3SJed Brown 1055f7a08781SBarry Smith static PetscErrorCode PCGASMSetSortIndices_GASM(PC pc,PetscBool doSort) 1056b1a0a8a3SJed Brown { 1057f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 1058b1a0a8a3SJed Brown 1059b1a0a8a3SJed Brown PetscFunctionBegin; 1060b1a0a8a3SJed Brown osm->sort_indices = doSort; 1061b1a0a8a3SJed Brown PetscFunctionReturn(0); 1062b1a0a8a3SJed Brown } 1063b1a0a8a3SJed Brown 106444fe18e5SDmitry Karpeev /* 10658f3b4b4dSDmitry Karpeev FIXME: This routine might need to be modified now that multiple ranks per subdomain are allowed. 106644fe18e5SDmitry Karpeev In particular, it would upset the global subdomain number calculation. 106744fe18e5SDmitry Karpeev */ 1068f7a08781SBarry Smith static PetscErrorCode PCGASMGetSubKSP_GASM(PC pc,PetscInt *n,PetscInt *first,KSP **ksp) 1069b1a0a8a3SJed Brown { 1070f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 1071b1a0a8a3SJed Brown 1072b1a0a8a3SJed Brown PetscFunctionBegin; 107308401ef6SPierre 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"); 1074b1a0a8a3SJed Brown 10752fa5cd67SKarl Rupp if (n) *n = osm->n; 1076ab3e923fSDmitry Karpeev if (first) { 10779566063dSJacob Faibussowitsch PetscCallMPI(MPI_Scan(&osm->n,first,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)pc))); 1078ab3e923fSDmitry Karpeev *first -= osm->n; 1079b1a0a8a3SJed Brown } 1080b1a0a8a3SJed Brown if (ksp) { 1081b1a0a8a3SJed Brown /* Assume that local solves are now different; not necessarily 108206b43e7eSDmitry Karpeev true, though! This flag is used only for PCView_GASM() */ 1083b1a0a8a3SJed Brown *ksp = osm->ksp; 10846a4f0f73SDmitry Karpeev osm->same_subdomain_solvers = PETSC_FALSE; 1085b1a0a8a3SJed Brown } 1086b1a0a8a3SJed Brown PetscFunctionReturn(0); 1087ab3e923fSDmitry Karpeev } /* PCGASMGetSubKSP_GASM() */ 1088b1a0a8a3SJed Brown 1089b1a0a8a3SJed Brown /*@C 109006b43e7eSDmitry Karpeev PCGASMSetSubdomains - Sets the subdomains for this processor 109106b43e7eSDmitry Karpeev for the additive Schwarz preconditioner. 1092b1a0a8a3SJed Brown 10938f3b4b4dSDmitry Karpeev Collective on pc 1094b1a0a8a3SJed Brown 1095b1a0a8a3SJed Brown Input Parameters: 10968f3b4b4dSDmitry Karpeev + pc - the preconditioner object 109706b43e7eSDmitry Karpeev . n - the number of subdomains for this processor 10988f3b4b4dSDmitry Karpeev . iis - the index sets that define the inner subdomains (or NULL for PETSc to determine subdomains) 10998f3b4b4dSDmitry 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) 1100b1a0a8a3SJed Brown 1101b1a0a8a3SJed Brown Notes: 110206b43e7eSDmitry Karpeev The IS indices use the parallel, global numbering of the vector entries. 11036a4f0f73SDmitry Karpeev Inner subdomains are those where the correction is applied. 11046a4f0f73SDmitry Karpeev Outer subdomains are those where the residual necessary to obtain the 11056a4f0f73SDmitry Karpeev corrections is obtained (see PCGASMType for the use of inner/outer subdomains). 11066a4f0f73SDmitry Karpeev Both inner and outer subdomains can extend over several processors. 11076a4f0f73SDmitry Karpeev This processor's portion of a subdomain is known as a local subdomain. 11086a4f0f73SDmitry Karpeev 1109670417b2SFande Kong Inner subdomains can not overlap with each other, do not have any entities from remote processors, 1110670417b2SFande Kong and have to cover the entire local subdomain owned by the current processor. The index sets on each 1111670417b2SFande Kong process should be ordered such that the ith local subdomain is connected to the ith remote subdomain 1112670417b2SFande Kong on another MPI process. 1113670417b2SFande Kong 11146a4f0f73SDmitry Karpeev By default the GASM preconditioner uses 1 (local) subdomain per processor. 11156a4f0f73SDmitry Karpeev 1116b1a0a8a3SJed Brown Level: advanced 1117b1a0a8a3SJed Brown 1118db781477SPatrick Sanan .seealso: `PCGASMSetOverlap()`, `PCGASMGetSubKSP()`, 1119db781477SPatrick Sanan `PCGASMCreateSubdomains2D()`, `PCGASMGetSubdomains()` 1120b1a0a8a3SJed Brown @*/ 11216a4f0f73SDmitry Karpeev PetscErrorCode PCGASMSetSubdomains(PC pc,PetscInt n,IS iis[],IS ois[]) 1122b1a0a8a3SJed Brown { 11238f3b4b4dSDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 1124b1a0a8a3SJed Brown 1125b1a0a8a3SJed Brown PetscFunctionBegin; 1126b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1127cac4c232SBarry Smith PetscTryMethod(pc,"PCGASMSetSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,n,iis,ois)); 11288f3b4b4dSDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 1129b1a0a8a3SJed Brown PetscFunctionReturn(0); 1130b1a0a8a3SJed Brown } 1131b1a0a8a3SJed Brown 1132b1a0a8a3SJed Brown /*@ 1133f746d493SDmitry Karpeev PCGASMSetOverlap - Sets the overlap between a pair of subdomains for the 1134b1a0a8a3SJed Brown additive Schwarz preconditioner. Either all or no processors in the 11358f3b4b4dSDmitry Karpeev pc communicator must call this routine. 1136b1a0a8a3SJed Brown 11378f3b4b4dSDmitry Karpeev Logically Collective on pc 1138b1a0a8a3SJed Brown 1139b1a0a8a3SJed Brown Input Parameters: 1140b1a0a8a3SJed Brown + pc - the preconditioner context 11418f3b4b4dSDmitry Karpeev - ovl - the amount of overlap between subdomains (ovl >= 0, default value = 0) 1142b1a0a8a3SJed Brown 1143b1a0a8a3SJed Brown Options Database Key: 114406b43e7eSDmitry Karpeev . -pc_gasm_overlap <overlap> - Sets overlap 1145b1a0a8a3SJed Brown 1146b1a0a8a3SJed Brown Notes: 114706b43e7eSDmitry Karpeev By default the GASM preconditioner uses 1 subdomain per processor. To use 11488f3b4b4dSDmitry Karpeev multiple subdomain per perocessor or "straddling" subdomains that intersect 11498f3b4b4dSDmitry Karpeev multiple processors use PCGASMSetSubdomains() (or option -pc_gasm_total_subdomains <n>). 1150b1a0a8a3SJed Brown 11518f3b4b4dSDmitry Karpeev The overlap defaults to 0, so if one desires that no additional 1152b1a0a8a3SJed Brown overlap be computed beyond what may have been set with a call to 11538f3b4b4dSDmitry Karpeev PCGASMSetSubdomains(), then ovl must be set to be 0. In particular, if one does 11548f3b4b4dSDmitry Karpeev not explicitly set the subdomains in application code, then all overlap would be computed 1155f746d493SDmitry Karpeev internally by PETSc, and using an overlap of 0 would result in an GASM 1156b1a0a8a3SJed Brown variant that is equivalent to the block Jacobi preconditioner. 1157b1a0a8a3SJed Brown 1158b1a0a8a3SJed Brown Note that one can define initial index sets with any overlap via 115906b43e7eSDmitry Karpeev PCGASMSetSubdomains(); the routine PCGASMSetOverlap() merely allows 116006b43e7eSDmitry Karpeev PETSc to extend that overlap further, if desired. 1161b1a0a8a3SJed Brown 1162b1a0a8a3SJed Brown Level: intermediate 1163b1a0a8a3SJed Brown 1164db781477SPatrick Sanan .seealso: `PCGASMSetSubdomains()`, `PCGASMGetSubKSP()`, 1165db781477SPatrick Sanan `PCGASMCreateSubdomains2D()`, `PCGASMGetSubdomains()` 1166b1a0a8a3SJed Brown @*/ 11677087cfbeSBarry Smith PetscErrorCode PCGASMSetOverlap(PC pc,PetscInt ovl) 1168b1a0a8a3SJed Brown { 11698f3b4b4dSDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 1170b1a0a8a3SJed Brown 1171b1a0a8a3SJed Brown PetscFunctionBegin; 1172b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1173b1a0a8a3SJed Brown PetscValidLogicalCollectiveInt(pc,ovl,2); 1174cac4c232SBarry Smith PetscTryMethod(pc,"PCGASMSetOverlap_C",(PC,PetscInt),(pc,ovl)); 11758f3b4b4dSDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 1176b1a0a8a3SJed Brown PetscFunctionReturn(0); 1177b1a0a8a3SJed Brown } 1178b1a0a8a3SJed Brown 1179b1a0a8a3SJed Brown /*@ 1180f746d493SDmitry Karpeev PCGASMSetType - Sets the type of restriction and interpolation used 1181b1a0a8a3SJed Brown for local problems in the additive Schwarz method. 1182b1a0a8a3SJed Brown 1183b1a0a8a3SJed Brown Logically Collective on PC 1184b1a0a8a3SJed Brown 1185b1a0a8a3SJed Brown Input Parameters: 1186b1a0a8a3SJed Brown + pc - the preconditioner context 1187f746d493SDmitry Karpeev - type - variant of GASM, one of 1188b1a0a8a3SJed Brown .vb 1189f746d493SDmitry Karpeev PC_GASM_BASIC - full interpolation and restriction 1190f746d493SDmitry Karpeev PC_GASM_RESTRICT - full restriction, local processor interpolation 1191f746d493SDmitry Karpeev PC_GASM_INTERPOLATE - full interpolation, local processor restriction 1192f746d493SDmitry Karpeev PC_GASM_NONE - local processor restriction and interpolation 1193b1a0a8a3SJed Brown .ve 1194b1a0a8a3SJed Brown 1195b1a0a8a3SJed Brown Options Database Key: 1196f746d493SDmitry Karpeev . -pc_gasm_type [basic,restrict,interpolate,none] - Sets GASM type 1197b1a0a8a3SJed Brown 1198b1a0a8a3SJed Brown Level: intermediate 1199b1a0a8a3SJed Brown 1200db781477SPatrick Sanan .seealso: `PCGASMSetSubdomains()`, `PCGASMGetSubKSP()`, 1201db781477SPatrick Sanan `PCGASMCreateSubdomains2D()` 1202b1a0a8a3SJed Brown @*/ 12037087cfbeSBarry Smith PetscErrorCode PCGASMSetType(PC pc,PCGASMType type) 1204b1a0a8a3SJed Brown { 1205b1a0a8a3SJed Brown PetscFunctionBegin; 1206b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1207b1a0a8a3SJed Brown PetscValidLogicalCollectiveEnum(pc,type,2); 1208cac4c232SBarry Smith PetscTryMethod(pc,"PCGASMSetType_C",(PC,PCGASMType),(pc,type)); 1209b1a0a8a3SJed Brown PetscFunctionReturn(0); 1210b1a0a8a3SJed Brown } 1211b1a0a8a3SJed Brown 1212b1a0a8a3SJed Brown /*@ 1213f746d493SDmitry Karpeev PCGASMSetSortIndices - Determines whether subdomain indices are sorted. 1214b1a0a8a3SJed Brown 1215b1a0a8a3SJed Brown Logically Collective on PC 1216b1a0a8a3SJed Brown 1217b1a0a8a3SJed Brown Input Parameters: 1218b1a0a8a3SJed Brown + pc - the preconditioner context 1219b1a0a8a3SJed Brown - doSort - sort the subdomain indices 1220b1a0a8a3SJed Brown 1221b1a0a8a3SJed Brown Level: intermediate 1222b1a0a8a3SJed Brown 1223db781477SPatrick Sanan .seealso: `PCGASMSetSubdomains()`, `PCGASMGetSubKSP()`, 1224db781477SPatrick Sanan `PCGASMCreateSubdomains2D()` 1225b1a0a8a3SJed Brown @*/ 12267087cfbeSBarry Smith PetscErrorCode PCGASMSetSortIndices(PC pc,PetscBool doSort) 1227b1a0a8a3SJed Brown { 1228b1a0a8a3SJed Brown PetscFunctionBegin; 1229b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1230b1a0a8a3SJed Brown PetscValidLogicalCollectiveBool(pc,doSort,2); 1231cac4c232SBarry Smith PetscTryMethod(pc,"PCGASMSetSortIndices_C",(PC,PetscBool),(pc,doSort)); 1232b1a0a8a3SJed Brown PetscFunctionReturn(0); 1233b1a0a8a3SJed Brown } 1234b1a0a8a3SJed Brown 1235b1a0a8a3SJed Brown /*@C 1236f746d493SDmitry Karpeev PCGASMGetSubKSP - Gets the local KSP contexts for all blocks on 1237b1a0a8a3SJed Brown this processor. 1238b1a0a8a3SJed Brown 1239b1a0a8a3SJed Brown Collective on PC iff first_local is requested 1240b1a0a8a3SJed Brown 1241b1a0a8a3SJed Brown Input Parameter: 1242b1a0a8a3SJed Brown . pc - the preconditioner context 1243b1a0a8a3SJed Brown 1244b1a0a8a3SJed Brown Output Parameters: 12450298fd71SBarry Smith + n_local - the number of blocks on this processor or NULL 12460298fd71SBarry Smith . first_local - the global number of the first block on this processor or NULL, 12470298fd71SBarry Smith all processors must request or all must pass NULL 1248b1a0a8a3SJed Brown - ksp - the array of KSP contexts 1249b1a0a8a3SJed Brown 1250b1a0a8a3SJed Brown Note: 1251f746d493SDmitry Karpeev After PCGASMGetSubKSP() the array of KSPes is not to be freed 1252b1a0a8a3SJed Brown 1253b1a0a8a3SJed Brown Currently for some matrix implementations only 1 block per processor 1254b1a0a8a3SJed Brown is supported. 1255b1a0a8a3SJed Brown 1256f746d493SDmitry Karpeev You must call KSPSetUp() before calling PCGASMGetSubKSP(). 1257b1a0a8a3SJed Brown 1258b1a0a8a3SJed Brown Level: advanced 1259b1a0a8a3SJed Brown 1260db781477SPatrick Sanan .seealso: `PCGASMSetSubdomains()`, `PCGASMSetOverlap()`, 1261db781477SPatrick Sanan `PCGASMCreateSubdomains2D()`, 1262b1a0a8a3SJed Brown @*/ 12637087cfbeSBarry Smith PetscErrorCode PCGASMGetSubKSP(PC pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp[]) 1264b1a0a8a3SJed Brown { 1265b1a0a8a3SJed Brown PetscFunctionBegin; 1266b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1267cac4c232SBarry Smith PetscUseMethod(pc,"PCGASMGetSubKSP_C",(PC,PetscInt*,PetscInt*,KSP **),(pc,n_local,first_local,ksp)); 1268b1a0a8a3SJed Brown PetscFunctionReturn(0); 1269b1a0a8a3SJed Brown } 1270b1a0a8a3SJed Brown 1271b1a0a8a3SJed Brown /* -------------------------------------------------------------------------------------*/ 1272b1a0a8a3SJed Brown /*MC 1273f746d493SDmitry Karpeev PCGASM - Use the (restricted) additive Schwarz method, each block is (approximately) solved with 1274b1a0a8a3SJed Brown its own KSP object. 1275b1a0a8a3SJed Brown 1276b1a0a8a3SJed Brown Options Database Keys: 12778f3b4b4dSDmitry Karpeev + -pc_gasm_total_subdomains <n> - Sets total number of local subdomains to be distributed among processors 127806b43e7eSDmitry Karpeev . -pc_gasm_view_subdomains - activates the printing of subdomain indices in PCView(), -ksp_view or -snes_view 127906b43e7eSDmitry Karpeev . -pc_gasm_print_subdomains - activates the printing of subdomain indices in PCSetUp() 128006b43e7eSDmitry Karpeev . -pc_gasm_overlap <ovl> - Sets overlap by which to (automatically) extend local subdomains 1281f746d493SDmitry Karpeev - -pc_gasm_type [basic,restrict,interpolate,none] - Sets GASM type 1282b1a0a8a3SJed Brown 1283b1a0a8a3SJed Brown IMPORTANT: If you run with, for example, 3 blocks on 1 processor or 3 blocks on 3 processors you 1284f746d493SDmitry Karpeev will get a different convergence rate due to the default option of -pc_gasm_type restrict. Use 1285f746d493SDmitry Karpeev -pc_gasm_type basic to use the standard GASM. 1286b1a0a8a3SJed Brown 128795452b02SPatrick Sanan Notes: 128895452b02SPatrick Sanan Blocks can be shared by multiple processes. 1289b1a0a8a3SJed Brown 1290b1a0a8a3SJed Brown To set options on the solvers for each block append -sub_ to all the KSP, and PC 1291b1a0a8a3SJed Brown options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly 1292b1a0a8a3SJed Brown 1293f746d493SDmitry Karpeev To set the options on the solvers separate for each block call PCGASMGetSubKSP() 1294b1a0a8a3SJed Brown and set the options directly on the resulting KSP object (you can access its PC 1295b1a0a8a3SJed Brown with KSPGetPC()) 1296b1a0a8a3SJed Brown 1297b1a0a8a3SJed Brown Level: beginner 1298b1a0a8a3SJed Brown 1299b1a0a8a3SJed Brown References: 1300606c0280SSatish Balay + * - M Dryja, OB Widlund, An additive variant of the Schwarz alternating method for the case of many subregions 130196a0c994SBarry Smith Courant Institute, New York University Technical report 1302606c0280SSatish Balay - * - Barry Smith, Petter Bjorstad, and William Gropp, Domain Decompositions: Parallel Multilevel Methods for Elliptic Partial Differential Equations, 130396a0c994SBarry Smith Cambridge University Press. 1304b1a0a8a3SJed Brown 1305db781477SPatrick Sanan .seealso: `PCCreate()`, `PCSetType()`, `PCType`, `PC`, 1306db781477SPatrick Sanan `PCBJACOBI`, `PCGASMGetSubKSP()`, `PCGASMSetSubdomains()`, 1307db781477SPatrick Sanan `PCSetModifySubMatrices()`, `PCGASMSetOverlap()`, `PCGASMSetType()` 1308b1a0a8a3SJed Brown 1309b1a0a8a3SJed Brown M*/ 1310b1a0a8a3SJed Brown 13118cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_GASM(PC pc) 1312b1a0a8a3SJed Brown { 1313f746d493SDmitry Karpeev PC_GASM *osm; 1314b1a0a8a3SJed Brown 1315b1a0a8a3SJed Brown PetscFunctionBegin; 13169566063dSJacob Faibussowitsch PetscCall(PetscNewLog(pc,&osm)); 13172fa5cd67SKarl Rupp 13188f3b4b4dSDmitry Karpeev osm->N = PETSC_DETERMINE; 131906b43e7eSDmitry Karpeev osm->n = PETSC_DECIDE; 13208f3b4b4dSDmitry Karpeev osm->nmax = PETSC_DETERMINE; 13218f3b4b4dSDmitry Karpeev osm->overlap = 0; 13220a545947SLisandro Dalcin osm->ksp = NULL; 13230a545947SLisandro Dalcin osm->gorestriction = NULL; 13240a545947SLisandro Dalcin osm->girestriction = NULL; 13250a545947SLisandro Dalcin osm->pctoouter = NULL; 13260a545947SLisandro Dalcin osm->gx = NULL; 13270a545947SLisandro Dalcin osm->gy = NULL; 13280a545947SLisandro Dalcin osm->x = NULL; 13290a545947SLisandro Dalcin osm->y = NULL; 13300a545947SLisandro Dalcin osm->pcx = NULL; 13310a545947SLisandro Dalcin osm->pcy = NULL; 13320a545947SLisandro Dalcin osm->permutationIS = NULL; 13330a545947SLisandro Dalcin osm->permutationP = NULL; 13340a545947SLisandro Dalcin osm->pcmat = NULL; 13350a545947SLisandro Dalcin osm->ois = NULL; 13360a545947SLisandro Dalcin osm->iis = NULL; 13370a545947SLisandro Dalcin osm->pmat = NULL; 1338f746d493SDmitry Karpeev osm->type = PC_GASM_RESTRICT; 13396a4f0f73SDmitry Karpeev osm->same_subdomain_solvers = PETSC_TRUE; 1340b1a0a8a3SJed Brown osm->sort_indices = PETSC_TRUE; 1341d709ea83SDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 1342ea91fabdSFande Kong osm->hierarchicalpartitioning = PETSC_FALSE; 1343b1a0a8a3SJed Brown 1344b1a0a8a3SJed Brown pc->data = (void*)osm; 1345f746d493SDmitry Karpeev pc->ops->apply = PCApply_GASM; 134648e38a8aSPierre Jolivet pc->ops->matapply = PCMatApply_GASM; 1347f746d493SDmitry Karpeev pc->ops->applytranspose = PCApplyTranspose_GASM; 1348f746d493SDmitry Karpeev pc->ops->setup = PCSetUp_GASM; 1349a06653b4SBarry Smith pc->ops->reset = PCReset_GASM; 1350f746d493SDmitry Karpeev pc->ops->destroy = PCDestroy_GASM; 1351f746d493SDmitry Karpeev pc->ops->setfromoptions = PCSetFromOptions_GASM; 1352f746d493SDmitry Karpeev pc->ops->setuponblocks = PCSetUpOnBlocks_GASM; 1353f746d493SDmitry Karpeev pc->ops->view = PCView_GASM; 13540a545947SLisandro Dalcin pc->ops->applyrichardson = NULL; 1355b1a0a8a3SJed Brown 13569566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetSubdomains_C",PCGASMSetSubdomains_GASM)); 13579566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetOverlap_C",PCGASMSetOverlap_GASM)); 13589566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetType_C",PCGASMSetType_GASM)); 13599566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetSortIndices_C",PCGASMSetSortIndices_GASM)); 13609566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)pc,"PCGASMGetSubKSP_C",PCGASMGetSubKSP_GASM)); 1361b1a0a8a3SJed Brown PetscFunctionReturn(0); 1362b1a0a8a3SJed Brown } 1363b1a0a8a3SJed Brown 13648f3b4b4dSDmitry Karpeev PetscErrorCode PCGASMCreateLocalSubdomains(Mat A, PetscInt nloc, IS *iis[]) 1365b1a0a8a3SJed Brown { 1366b1a0a8a3SJed Brown MatPartitioning mpart; 1367b1a0a8a3SJed Brown const char *prefix; 1368b1a0a8a3SJed Brown PetscInt i,j,rstart,rend,bs; 1369976e8c5aSLisandro Dalcin PetscBool hasop, isbaij = PETSC_FALSE,foundpart = PETSC_FALSE; 13700298fd71SBarry Smith Mat Ad = NULL, adj; 1371b1a0a8a3SJed Brown IS ispart,isnumb,*is; 1372b1a0a8a3SJed Brown 1373b1a0a8a3SJed Brown PetscFunctionBegin; 137463a3b9bcSJacob Faibussowitsch PetscCheck(nloc >= 1,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of local subdomains must > 0, got nloc = %" PetscInt_FMT,nloc); 1375b1a0a8a3SJed Brown 1376b1a0a8a3SJed Brown /* Get prefix, row distribution, and block size */ 13779566063dSJacob Faibussowitsch PetscCall(MatGetOptionsPrefix(A,&prefix)); 13789566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(A,&rstart,&rend)); 13799566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(A,&bs)); 13802472a847SBarry 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); 1381b1a0a8a3SJed Brown 1382b1a0a8a3SJed Brown /* Get diagonal block from matrix if possible */ 13839566063dSJacob Faibussowitsch PetscCall(MatHasOperation(A,MATOP_GET_DIAGONAL_BLOCK,&hasop)); 1384976e8c5aSLisandro Dalcin if (hasop) { 13859566063dSJacob Faibussowitsch PetscCall(MatGetDiagonalBlock(A,&Ad)); 1386b1a0a8a3SJed Brown } 1387b1a0a8a3SJed Brown if (Ad) { 13889566063dSJacob Faibussowitsch PetscCall(PetscObjectBaseTypeCompare((PetscObject)Ad,MATSEQBAIJ,&isbaij)); 13899566063dSJacob Faibussowitsch if (!isbaij) PetscCall(PetscObjectBaseTypeCompare((PetscObject)Ad,MATSEQSBAIJ,&isbaij)); 1390b1a0a8a3SJed Brown } 13918f3b4b4dSDmitry Karpeev if (Ad && nloc > 1) { 1392b1a0a8a3SJed Brown PetscBool match,done; 1393b1a0a8a3SJed Brown /* Try to setup a good matrix partitioning if available */ 13949566063dSJacob Faibussowitsch PetscCall(MatPartitioningCreate(PETSC_COMM_SELF,&mpart)); 13959566063dSJacob Faibussowitsch PetscCall(PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix)); 13969566063dSJacob Faibussowitsch PetscCall(MatPartitioningSetFromOptions(mpart)); 13979566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGCURRENT,&match)); 1398b1a0a8a3SJed Brown if (!match) { 13999566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGSQUARE,&match)); 1400b1a0a8a3SJed Brown } 1401b1a0a8a3SJed Brown if (!match) { /* assume a "good" partitioner is available */ 14021a83f524SJed Brown PetscInt na; 14031a83f524SJed Brown const PetscInt *ia,*ja; 14049566063dSJacob Faibussowitsch PetscCall(MatGetRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done)); 1405b1a0a8a3SJed Brown if (done) { 1406b1a0a8a3SJed Brown /* Build adjacency matrix by hand. Unfortunately a call to 1407b1a0a8a3SJed Brown MatConvert(Ad,MATMPIADJ,MAT_INITIAL_MATRIX,&adj) will 1408b1a0a8a3SJed Brown remove the block-aij structure and we cannot expect 1409b1a0a8a3SJed Brown MatPartitioning to split vertices as we need */ 14100a545947SLisandro Dalcin PetscInt i,j,len,nnz,cnt,*iia=NULL,*jja=NULL; 14111a83f524SJed Brown const PetscInt *row; 1412b1a0a8a3SJed Brown nnz = 0; 1413b1a0a8a3SJed Brown for (i=0; i<na; i++) { /* count number of nonzeros */ 1414b1a0a8a3SJed Brown len = ia[i+1] - ia[i]; 1415b1a0a8a3SJed Brown row = ja + ia[i]; 1416b1a0a8a3SJed Brown for (j=0; j<len; j++) { 1417b1a0a8a3SJed Brown if (row[j] == i) { /* don't count diagonal */ 1418b1a0a8a3SJed Brown len--; break; 1419b1a0a8a3SJed Brown } 1420b1a0a8a3SJed Brown } 1421b1a0a8a3SJed Brown nnz += len; 1422b1a0a8a3SJed Brown } 14239566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(na+1,&iia)); 14249566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nnz,&jja)); 1425b1a0a8a3SJed Brown nnz = 0; 1426b1a0a8a3SJed Brown iia[0] = 0; 1427b1a0a8a3SJed Brown for (i=0; i<na; i++) { /* fill adjacency */ 1428b1a0a8a3SJed Brown cnt = 0; 1429b1a0a8a3SJed Brown len = ia[i+1] - ia[i]; 1430b1a0a8a3SJed Brown row = ja + ia[i]; 1431b1a0a8a3SJed Brown for (j=0; j<len; j++) { 14322fa5cd67SKarl Rupp if (row[j] != i) jja[nnz+cnt++] = row[j]; /* if not diagonal */ 1433b1a0a8a3SJed Brown } 1434b1a0a8a3SJed Brown nnz += cnt; 1435b1a0a8a3SJed Brown iia[i+1] = nnz; 1436b1a0a8a3SJed Brown } 1437b1a0a8a3SJed Brown /* Partitioning of the adjacency matrix */ 14389566063dSJacob Faibussowitsch PetscCall(MatCreateMPIAdj(PETSC_COMM_SELF,na,na,iia,jja,NULL,&adj)); 14399566063dSJacob Faibussowitsch PetscCall(MatPartitioningSetAdjacency(mpart,adj)); 14409566063dSJacob Faibussowitsch PetscCall(MatPartitioningSetNParts(mpart,nloc)); 14419566063dSJacob Faibussowitsch PetscCall(MatPartitioningApply(mpart,&ispart)); 14429566063dSJacob Faibussowitsch PetscCall(ISPartitioningToNumbering(ispart,&isnumb)); 14439566063dSJacob Faibussowitsch PetscCall(MatDestroy(&adj)); 1444b1a0a8a3SJed Brown foundpart = PETSC_TRUE; 1445b1a0a8a3SJed Brown } 14469566063dSJacob Faibussowitsch PetscCall(MatRestoreRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done)); 1447b1a0a8a3SJed Brown } 14489566063dSJacob Faibussowitsch PetscCall(MatPartitioningDestroy(&mpart)); 1449b1a0a8a3SJed Brown } 14509566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nloc,&is)); 1451b1a0a8a3SJed Brown if (!foundpart) { 1452b1a0a8a3SJed Brown 1453b1a0a8a3SJed Brown /* Partitioning by contiguous chunks of rows */ 1454b1a0a8a3SJed Brown 1455b1a0a8a3SJed Brown PetscInt mbs = (rend-rstart)/bs; 1456b1a0a8a3SJed Brown PetscInt start = rstart; 14578f3b4b4dSDmitry Karpeev for (i=0; i<nloc; i++) { 14588f3b4b4dSDmitry Karpeev PetscInt count = (mbs/nloc + ((mbs % nloc) > i)) * bs; 14599566063dSJacob Faibussowitsch PetscCall(ISCreateStride(PETSC_COMM_SELF,count,start,1,&is[i])); 1460b1a0a8a3SJed Brown start += count; 1461b1a0a8a3SJed Brown } 1462b1a0a8a3SJed Brown 1463b1a0a8a3SJed Brown } else { 1464b1a0a8a3SJed Brown 1465b1a0a8a3SJed Brown /* Partitioning by adjacency of diagonal block */ 1466b1a0a8a3SJed Brown 1467b1a0a8a3SJed Brown const PetscInt *numbering; 1468b1a0a8a3SJed Brown PetscInt *count,nidx,*indices,*newidx,start=0; 1469b1a0a8a3SJed Brown /* Get node count in each partition */ 14709566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nloc,&count)); 14719566063dSJacob Faibussowitsch PetscCall(ISPartitioningCount(ispart,nloc,count)); 1472b1a0a8a3SJed Brown if (isbaij && bs > 1) { /* adjust for the block-aij case */ 14738f3b4b4dSDmitry Karpeev for (i=0; i<nloc; i++) count[i] *= bs; 1474b1a0a8a3SJed Brown } 1475b1a0a8a3SJed Brown /* Build indices from node numbering */ 14769566063dSJacob Faibussowitsch PetscCall(ISGetLocalSize(isnumb,&nidx)); 14779566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nidx,&indices)); 1478b1a0a8a3SJed Brown for (i=0; i<nidx; i++) indices[i] = i; /* needs to be initialized */ 14799566063dSJacob Faibussowitsch PetscCall(ISGetIndices(isnumb,&numbering)); 14809566063dSJacob Faibussowitsch PetscCall(PetscSortIntWithPermutation(nidx,numbering,indices)); 14819566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(isnumb,&numbering)); 1482b1a0a8a3SJed Brown if (isbaij && bs > 1) { /* adjust for the block-aij case */ 14839566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nidx*bs,&newidx)); 14842fa5cd67SKarl Rupp for (i=0; i<nidx; i++) { 14852fa5cd67SKarl Rupp for (j=0; j<bs; j++) newidx[i*bs+j] = indices[i]*bs + j; 14862fa5cd67SKarl Rupp } 14879566063dSJacob Faibussowitsch PetscCall(PetscFree(indices)); 1488b1a0a8a3SJed Brown nidx *= bs; 1489b1a0a8a3SJed Brown indices = newidx; 1490b1a0a8a3SJed Brown } 1491b1a0a8a3SJed Brown /* Shift to get global indices */ 1492b1a0a8a3SJed Brown for (i=0; i<nidx; i++) indices[i] += rstart; 1493b1a0a8a3SJed Brown 1494b1a0a8a3SJed Brown /* Build the index sets for each block */ 14958f3b4b4dSDmitry Karpeev for (i=0; i<nloc; i++) { 14969566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(PETSC_COMM_SELF,count[i],&indices[start],PETSC_COPY_VALUES,&is[i])); 14979566063dSJacob Faibussowitsch PetscCall(ISSort(is[i])); 1498b1a0a8a3SJed Brown start += count[i]; 1499b1a0a8a3SJed Brown } 1500b1a0a8a3SJed Brown 15019566063dSJacob Faibussowitsch PetscCall(PetscFree(count)); 15029566063dSJacob Faibussowitsch PetscCall(PetscFree(indices)); 15039566063dSJacob Faibussowitsch PetscCall(ISDestroy(&isnumb)); 15049566063dSJacob Faibussowitsch PetscCall(ISDestroy(&ispart)); 1505b1a0a8a3SJed Brown } 15066a4f0f73SDmitry Karpeev *iis = is; 15078f3b4b4dSDmitry Karpeev PetscFunctionReturn(0); 15088f3b4b4dSDmitry Karpeev } 15098f3b4b4dSDmitry Karpeev 1510b541e6a4SDmitry Karpeev PETSC_INTERN PetscErrorCode PCGASMCreateStraddlingSubdomains(Mat A,PetscInt N,PetscInt *n,IS *iis[]) 15118f3b4b4dSDmitry Karpeev { 15128f3b4b4dSDmitry Karpeev PetscFunctionBegin; 15139566063dSJacob Faibussowitsch PetscCall(MatSubdomainsCreateCoalesce(A,N,n,iis)); 15148f3b4b4dSDmitry Karpeev PetscFunctionReturn(0); 15158f3b4b4dSDmitry Karpeev } 15168f3b4b4dSDmitry Karpeev 15178f3b4b4dSDmitry Karpeev /*@C 15188f3b4b4dSDmitry Karpeev PCGASMCreateSubdomains - Creates n index sets defining n nonoverlapping subdomains for the additive 15198f3b4b4dSDmitry Karpeev Schwarz preconditioner for a any problem based on its matrix. 15208f3b4b4dSDmitry Karpeev 15218f3b4b4dSDmitry Karpeev Collective 15228f3b4b4dSDmitry Karpeev 15238f3b4b4dSDmitry Karpeev Input Parameters: 15248f3b4b4dSDmitry Karpeev + A - The global matrix operator 15258f3b4b4dSDmitry Karpeev - N - the number of global subdomains requested 15268f3b4b4dSDmitry Karpeev 15278f3b4b4dSDmitry Karpeev Output Parameters: 15288f3b4b4dSDmitry Karpeev + n - the number of subdomains created on this processor 15298f3b4b4dSDmitry Karpeev - iis - the array of index sets defining the local inner subdomains (on which the correction is applied) 15308f3b4b4dSDmitry Karpeev 15318f3b4b4dSDmitry Karpeev Level: advanced 15328f3b4b4dSDmitry Karpeev 15338f3b4b4dSDmitry Karpeev Note: When N >= A's communicator size, each subdomain is local -- contained within a single processor. 15348f3b4b4dSDmitry Karpeev When N < size, the subdomains are 'straddling' (processor boundaries) and are no longer local. 15358f3b4b4dSDmitry Karpeev The resulting subdomains can be use in PCGASMSetSubdomains(pc,n,iss,NULL). The overlapping 15368f3b4b4dSDmitry Karpeev outer subdomains will be automatically generated from these according to the requested amount of 15378f3b4b4dSDmitry Karpeev overlap; this is currently supported only with local subdomains. 15388f3b4b4dSDmitry Karpeev 1539db781477SPatrick Sanan .seealso: `PCGASMSetSubdomains()`, `PCGASMDestroySubdomains()` 15408f3b4b4dSDmitry Karpeev @*/ 1541b541e6a4SDmitry Karpeev PetscErrorCode PCGASMCreateSubdomains(Mat A,PetscInt N,PetscInt *n,IS *iis[]) 15428f3b4b4dSDmitry Karpeev { 15438f3b4b4dSDmitry Karpeev PetscMPIInt size; 15448f3b4b4dSDmitry Karpeev 15458f3b4b4dSDmitry Karpeev PetscFunctionBegin; 15468f3b4b4dSDmitry Karpeev PetscValidHeaderSpecific(A,MAT_CLASSID,1); 15478f3b4b4dSDmitry Karpeev PetscValidPointer(iis,4); 15488f3b4b4dSDmitry Karpeev 154963a3b9bcSJacob Faibussowitsch PetscCheck(N >= 1,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Number of subdomains must be > 0, N = %" PetscInt_FMT,N); 15509566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A),&size)); 15518f3b4b4dSDmitry Karpeev if (N >= size) { 15528f3b4b4dSDmitry Karpeev *n = N/size + (N%size); 15539566063dSJacob Faibussowitsch PetscCall(PCGASMCreateLocalSubdomains(A,*n,iis)); 15546a4f0f73SDmitry Karpeev } else { 15559566063dSJacob Faibussowitsch PetscCall(PCGASMCreateStraddlingSubdomains(A,N,n,iis)); 15566a4f0f73SDmitry Karpeev } 1557b1a0a8a3SJed Brown PetscFunctionReturn(0); 1558b1a0a8a3SJed Brown } 1559b1a0a8a3SJed Brown 1560b1a0a8a3SJed Brown /*@C 1561f746d493SDmitry Karpeev PCGASMDestroySubdomains - Destroys the index sets created with 15628f3b4b4dSDmitry Karpeev PCGASMCreateSubdomains() or PCGASMCreateSubdomains2D. Should be 156306b43e7eSDmitry Karpeev called after setting subdomains with PCGASMSetSubdomains(). 1564b1a0a8a3SJed Brown 1565b1a0a8a3SJed Brown Collective 1566b1a0a8a3SJed Brown 1567b1a0a8a3SJed Brown Input Parameters: 1568b1a0a8a3SJed Brown + n - the number of index sets 15696a4f0f73SDmitry Karpeev . iis - the array of inner subdomains, 15700298fd71SBarry Smith - ois - the array of outer subdomains, can be NULL 1571b1a0a8a3SJed Brown 15726a4f0f73SDmitry Karpeev Level: intermediate 15736a4f0f73SDmitry Karpeev 157495452b02SPatrick Sanan Notes: 157595452b02SPatrick Sanan this is merely a convenience subroutine that walks each list, 15762c112581SDmitry Karpeev destroys each IS on the list, and then frees the list. At the end the 15772c112581SDmitry Karpeev list pointers are set to NULL. 1578b1a0a8a3SJed Brown 1579db781477SPatrick Sanan .seealso: `PCGASMCreateSubdomains()`, `PCGASMSetSubdomains()` 1580b1a0a8a3SJed Brown @*/ 15812c112581SDmitry Karpeev PetscErrorCode PCGASMDestroySubdomains(PetscInt n,IS **iis,IS **ois) 1582b1a0a8a3SJed Brown { 1583b1a0a8a3SJed Brown PetscInt i; 15845fd66863SKarl Rupp 1585b1a0a8a3SJed Brown PetscFunctionBegin; 1586a35b7b57SDmitry Karpeev if (n <= 0) PetscFunctionReturn(0); 15876a4f0f73SDmitry Karpeev if (ois) { 15882c112581SDmitry Karpeev PetscValidPointer(ois,3); 15892c112581SDmitry Karpeev if (*ois) { 15902c112581SDmitry Karpeev PetscValidPointer(*ois,3); 1591a35b7b57SDmitry Karpeev for (i=0; i<n; i++) { 15929566063dSJacob Faibussowitsch PetscCall(ISDestroy(&(*ois)[i])); 1593a35b7b57SDmitry Karpeev } 15949566063dSJacob Faibussowitsch PetscCall(PetscFree((*ois))); 15952c112581SDmitry Karpeev } 1596b1a0a8a3SJed Brown } 1597b9d0fdaaSFande Kong if (iis) { 1598b9d0fdaaSFande Kong PetscValidPointer(iis,2); 1599b9d0fdaaSFande Kong if (*iis) { 1600b9d0fdaaSFande Kong PetscValidPointer(*iis,2); 1601b9d0fdaaSFande Kong for (i=0; i<n; i++) { 16029566063dSJacob Faibussowitsch PetscCall(ISDestroy(&(*iis)[i])); 1603b9d0fdaaSFande Kong } 16049566063dSJacob Faibussowitsch PetscCall(PetscFree((*iis))); 1605b9d0fdaaSFande Kong } 1606b9d0fdaaSFande Kong } 1607b1a0a8a3SJed Brown PetscFunctionReturn(0); 1608b1a0a8a3SJed Brown } 1609b1a0a8a3SJed Brown 1610af538404SDmitry Karpeev #define PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,xleft_loc,ylow_loc,xright_loc,yhigh_loc,n) \ 1611af538404SDmitry Karpeev { \ 1612af538404SDmitry Karpeev PetscInt first_row = first/M, last_row = last/M+1; \ 1613af538404SDmitry Karpeev /* \ 1614af538404SDmitry Karpeev Compute ylow_loc and yhigh_loc so that (ylow_loc,xleft) and (yhigh_loc,xright) are the corners \ 1615af538404SDmitry Karpeev of the bounding box of the intersection of the subdomain with the local ownership range (local \ 1616af538404SDmitry Karpeev subdomain). \ 1617af538404SDmitry Karpeev Also compute xleft_loc and xright_loc as the lower and upper bounds on the first and last rows \ 1618af538404SDmitry Karpeev of the intersection. \ 1619af538404SDmitry Karpeev */ \ 1620af538404SDmitry Karpeev /* ylow_loc is the grid row containing the first element of the local sumbdomain */ \ 1621eec7e3faSDmitry Karpeev *ylow_loc = PetscMax(first_row,ylow); \ 1622af538404SDmitry Karpeev /* xleft_loc is the offset of first element of the local subdomain within its grid row (might actually be outside the local subdomain) */ \ 1623eec7e3faSDmitry Karpeev *xleft_loc = *ylow_loc==first_row ? PetscMax(first%M,xleft) : xleft; \ 1624af538404SDmitry Karpeev /* yhigh_loc is the grid row above the last local subdomain element */ \ 1625eec7e3faSDmitry Karpeev *yhigh_loc = PetscMin(last_row,yhigh); \ 1626af538404SDmitry Karpeev /* xright is the offset of the end of the local subdomain within its grid row (might actually be outside the local subdomain) */ \ 1627eec7e3faSDmitry Karpeev *xright_loc = *yhigh_loc==last_row ? PetscMin(xright,last%M) : xright; \ 1628af538404SDmitry Karpeev /* Now compute the size of the local subdomain n. */ \ 1629c3518bceSDmitry Karpeev *n = 0; \ 1630eec7e3faSDmitry Karpeev if (*ylow_loc < *yhigh_loc) { \ 1631af538404SDmitry Karpeev PetscInt width = xright-xleft; \ 1632eec7e3faSDmitry Karpeev *n += width*(*yhigh_loc-*ylow_loc-1); \ 1633eec7e3faSDmitry Karpeev *n += PetscMin(PetscMax(*xright_loc-xleft,0),width); \ 1634eec7e3faSDmitry Karpeev *n -= PetscMin(PetscMax(*xleft_loc-xleft,0), width); \ 1635af538404SDmitry Karpeev } \ 1636af538404SDmitry Karpeev } 1637af538404SDmitry Karpeev 1638b1a0a8a3SJed Brown /*@ 1639f746d493SDmitry Karpeev PCGASMCreateSubdomains2D - Creates the index sets for the overlapping Schwarz 1640b1a0a8a3SJed Brown preconditioner for a two-dimensional problem on a regular grid. 1641b1a0a8a3SJed Brown 1642af538404SDmitry Karpeev Collective 1643b1a0a8a3SJed Brown 1644b1a0a8a3SJed Brown Input Parameters: 16456b867d5aSJose E. Roman + pc - the preconditioner context 16466b867d5aSJose E. Roman . M - the global number of grid points in the x direction 16476b867d5aSJose E. Roman . N - the global number of grid points in the y direction 16486b867d5aSJose E. Roman . Mdomains - the global number of subdomains in the x direction 16496b867d5aSJose E. Roman . Ndomains - the global number of subdomains in the y direction 1650b1a0a8a3SJed Brown . dof - degrees of freedom per node 1651b1a0a8a3SJed Brown - overlap - overlap in mesh lines 1652b1a0a8a3SJed Brown 1653b1a0a8a3SJed Brown Output Parameters: 1654af538404SDmitry Karpeev + Nsub - the number of local subdomains created 16556a4f0f73SDmitry Karpeev . iis - array of index sets defining inner (nonoverlapping) subdomains 16566a4f0f73SDmitry Karpeev - ois - array of index sets defining outer (overlapping, if overlap > 0) subdomains 1657b1a0a8a3SJed Brown 1658b1a0a8a3SJed Brown Level: advanced 1659b1a0a8a3SJed Brown 1660db781477SPatrick Sanan .seealso: `PCGASMSetSubdomains()`, `PCGASMGetSubKSP()`, `PCGASMSetOverlap()` 1661b1a0a8a3SJed Brown @*/ 16626a4f0f73SDmitry Karpeev PetscErrorCode PCGASMCreateSubdomains2D(PC pc,PetscInt M,PetscInt N,PetscInt Mdomains,PetscInt Ndomains,PetscInt dof,PetscInt overlap,PetscInt *nsub,IS **iis,IS **ois) 1663b1a0a8a3SJed Brown { 16642bbb417fSDmitry Karpeev PetscMPIInt size, rank; 16652bbb417fSDmitry Karpeev PetscInt i, j; 16662bbb417fSDmitry Karpeev PetscInt maxheight, maxwidth; 16672bbb417fSDmitry Karpeev PetscInt xstart, xleft, xright, xleft_loc, xright_loc; 16682bbb417fSDmitry Karpeev PetscInt ystart, ylow, yhigh, ylow_loc, yhigh_loc; 1669eec7e3faSDmitry Karpeev PetscInt x[2][2], y[2][2], n[2]; 16702bbb417fSDmitry Karpeev PetscInt first, last; 16712bbb417fSDmitry Karpeev PetscInt nidx, *idx; 16722bbb417fSDmitry Karpeev PetscInt ii,jj,s,q,d; 1673af538404SDmitry Karpeev PetscInt k,kk; 16742bbb417fSDmitry Karpeev PetscMPIInt color; 16752bbb417fSDmitry Karpeev MPI_Comm comm, subcomm; 16760a545947SLisandro Dalcin IS **xis = NULL, **is = ois, **is_local = iis; 1677d34fcf5fSBarry Smith 1678b1a0a8a3SJed Brown PetscFunctionBegin; 16799566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)pc, &comm)); 16809566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(comm, &size)); 16819566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank)); 16829566063dSJacob Faibussowitsch PetscCall(MatGetOwnershipRange(pc->pmat, &first, &last)); 16832472a847SBarry 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 ") " 168463a3b9bcSJacob Faibussowitsch "does not respect the number of degrees of freedom per grid point %" PetscInt_FMT, first, last, dof); 1685d34fcf5fSBarry Smith 1686af538404SDmitry Karpeev /* Determine the number of domains with nonzero intersections with the local ownership range. */ 16872bbb417fSDmitry Karpeev s = 0; 1688af538404SDmitry Karpeev ystart = 0; 1689af538404SDmitry Karpeev for (j=0; j<Ndomains; ++j) { 1690af538404SDmitry Karpeev maxheight = N/Ndomains + ((N % Ndomains) > j); /* Maximal height of subdomain */ 169163a3b9bcSJacob 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); 1692eec7e3faSDmitry Karpeev /* Vertical domain limits with an overlap. */ 1693eec7e3faSDmitry Karpeev ylow = PetscMax(ystart - overlap,0); 1694eec7e3faSDmitry Karpeev yhigh = PetscMin(ystart + maxheight + overlap,N); 1695b1a0a8a3SJed Brown xstart = 0; 1696af538404SDmitry Karpeev for (i=0; i<Mdomains; ++i) { 1697af538404SDmitry Karpeev maxwidth = M/Mdomains + ((M % Mdomains) > i); /* Maximal width of subdomain */ 169863a3b9bcSJacob 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); 1699eec7e3faSDmitry Karpeev /* Horizontal domain limits with an overlap. */ 1700eec7e3faSDmitry Karpeev xleft = PetscMax(xstart - overlap,0); 1701eec7e3faSDmitry Karpeev xright = PetscMin(xstart + maxwidth + overlap,M); 1702af538404SDmitry Karpeev /* 1703af538404SDmitry Karpeev Determine whether this subdomain intersects this processor's ownership range of pc->pmat. 1704af538404SDmitry Karpeev */ 1705c3518bceSDmitry Karpeev PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,(&xleft_loc),(&ylow_loc),(&xright_loc),(&yhigh_loc),(&nidx)); 17062fa5cd67SKarl Rupp if (nidx) ++s; 1707af538404SDmitry Karpeev xstart += maxwidth; 1708af538404SDmitry Karpeev } /* for (i = 0; i < Mdomains; ++i) */ 1709af538404SDmitry Karpeev ystart += maxheight; 1710af538404SDmitry Karpeev } /* for (j = 0; j < Ndomains; ++j) */ 17112fa5cd67SKarl Rupp 1712af538404SDmitry Karpeev /* Now we can allocate the necessary number of ISs. */ 1713af538404SDmitry Karpeev *nsub = s; 17149566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(*nsub,is)); 17159566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(*nsub,is_local)); 1716af538404SDmitry Karpeev s = 0; 1717af538404SDmitry Karpeev ystart = 0; 1718af538404SDmitry Karpeev for (j=0; j<Ndomains; ++j) { 1719af538404SDmitry Karpeev maxheight = N/Ndomains + ((N % Ndomains) > j); /* Maximal height of subdomain */ 172063a3b9bcSJacob 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); 1721af538404SDmitry Karpeev /* Vertical domain limits with an overlap. */ 1722eec7e3faSDmitry Karpeev y[0][0] = PetscMax(ystart - overlap,0); 1723eec7e3faSDmitry Karpeev y[0][1] = PetscMin(ystart + maxheight + overlap,N); 1724eec7e3faSDmitry Karpeev /* Vertical domain limits without an overlap. */ 1725eec7e3faSDmitry Karpeev y[1][0] = ystart; 1726eec7e3faSDmitry Karpeev y[1][1] = PetscMin(ystart + maxheight,N); 1727eec7e3faSDmitry Karpeev xstart = 0; 1728af538404SDmitry Karpeev for (i=0; i<Mdomains; ++i) { 1729af538404SDmitry Karpeev maxwidth = M/Mdomains + ((M % Mdomains) > i); /* Maximal width of subdomain */ 173063a3b9bcSJacob 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); 1731af538404SDmitry Karpeev /* Horizontal domain limits with an overlap. */ 1732eec7e3faSDmitry Karpeev x[0][0] = PetscMax(xstart - overlap,0); 1733eec7e3faSDmitry Karpeev x[0][1] = PetscMin(xstart + maxwidth + overlap,M); 1734eec7e3faSDmitry Karpeev /* Horizontal domain limits without an overlap. */ 1735eec7e3faSDmitry Karpeev x[1][0] = xstart; 1736eec7e3faSDmitry Karpeev x[1][1] = PetscMin(xstart+maxwidth,M); 17372bbb417fSDmitry Karpeev /* 17382bbb417fSDmitry Karpeev Determine whether this domain intersects this processor's ownership range of pc->pmat. 17392bbb417fSDmitry Karpeev Do this twice: first for the domains with overlaps, and once without. 17402bbb417fSDmitry Karpeev During the first pass create the subcommunicators, and use them on the second pass as well. 17412bbb417fSDmitry Karpeev */ 17422bbb417fSDmitry Karpeev for (q = 0; q < 2; ++q) { 1743cc96b2e5SBarry Smith PetscBool split = PETSC_FALSE; 17442bbb417fSDmitry Karpeev /* 17452bbb417fSDmitry Karpeev domain limits, (xleft, xright) and (ylow, yheigh) are adjusted 17462bbb417fSDmitry Karpeev according to whether the domain with an overlap or without is considered. 17472bbb417fSDmitry Karpeev */ 17482bbb417fSDmitry Karpeev xleft = x[q][0]; xright = x[q][1]; 17492bbb417fSDmitry Karpeev ylow = y[q][0]; yhigh = y[q][1]; 1750c3518bceSDmitry Karpeev PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,(&xleft_loc),(&ylow_loc),(&xright_loc),(&yhigh_loc),(&nidx)); 1751af538404SDmitry Karpeev nidx *= dof; 1752eec7e3faSDmitry Karpeev n[q] = nidx; 17532bbb417fSDmitry Karpeev /* 1754eec7e3faSDmitry Karpeev Based on the counted number of indices in the local domain *with an overlap*, 1755af538404SDmitry Karpeev construct a subcommunicator of all the processors supporting this domain. 1756eec7e3faSDmitry Karpeev Observe that a domain with an overlap might have nontrivial local support, 1757eec7e3faSDmitry Karpeev while the domain without an overlap might not. Hence, the decision to participate 1758eec7e3faSDmitry Karpeev in the subcommunicator must be based on the domain with an overlap. 17592bbb417fSDmitry Karpeev */ 17602bbb417fSDmitry Karpeev if (q == 0) { 17612fa5cd67SKarl Rupp if (nidx) color = 1; 17622fa5cd67SKarl Rupp else color = MPI_UNDEFINED; 17639566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_split(comm, color, rank, &subcomm)); 1764cc96b2e5SBarry Smith split = PETSC_TRUE; 17652bbb417fSDmitry Karpeev } 1766af538404SDmitry Karpeev /* 1767eec7e3faSDmitry Karpeev Proceed only if the number of local indices *with an overlap* is nonzero. 1768af538404SDmitry Karpeev */ 1769eec7e3faSDmitry Karpeev if (n[0]) { 17702fa5cd67SKarl Rupp if (q == 0) xis = is; 1771af538404SDmitry Karpeev if (q == 1) { 1772af538404SDmitry Karpeev /* 1773eec7e3faSDmitry Karpeev The IS for the no-overlap subdomain shares a communicator with the overlapping domain. 1774af538404SDmitry Karpeev Moreover, if the overlap is zero, the two ISs are identical. 1775af538404SDmitry Karpeev */ 1776b1a0a8a3SJed Brown if (overlap == 0) { 1777eec7e3faSDmitry Karpeev (*is_local)[s] = (*is)[s]; 17789566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)(*is)[s])); 17792bbb417fSDmitry Karpeev continue; 1780d34fcf5fSBarry Smith } else { 17816a4f0f73SDmitry Karpeev xis = is_local; 1782eec7e3faSDmitry Karpeev subcomm = ((PetscObject)(*is)[s])->comm; 17832bbb417fSDmitry Karpeev } 1784af538404SDmitry Karpeev } /* if (q == 1) */ 17850298fd71SBarry Smith idx = NULL; 17869566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(nidx,&idx)); 1787eec7e3faSDmitry Karpeev if (nidx) { 17882bbb417fSDmitry Karpeev k = 0; 17892bbb417fSDmitry Karpeev for (jj=ylow_loc; jj<yhigh_loc; ++jj) { 1790af538404SDmitry Karpeev PetscInt x0 = (jj==ylow_loc) ? xleft_loc : xleft; 1791af538404SDmitry Karpeev PetscInt x1 = (jj==yhigh_loc-1) ? xright_loc : xright; 1792af538404SDmitry Karpeev kk = dof*(M*jj + x0); 17932bbb417fSDmitry Karpeev for (ii=x0; ii<x1; ++ii) { 17942bbb417fSDmitry Karpeev for (d = 0; d < dof; ++d) { 17952bbb417fSDmitry Karpeev idx[k++] = kk++; 1796b1a0a8a3SJed Brown } 1797b1a0a8a3SJed Brown } 1798b1a0a8a3SJed Brown } 1799eec7e3faSDmitry Karpeev } 18009566063dSJacob Faibussowitsch PetscCall(ISCreateGeneral(subcomm,nidx,idx,PETSC_OWN_POINTER,(*xis)+s)); 1801cc96b2e5SBarry Smith if (split) { 18029566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_free(&subcomm)); 1803cc96b2e5SBarry Smith } 1804eec7e3faSDmitry Karpeev }/* if (n[0]) */ 18052bbb417fSDmitry Karpeev }/* for (q = 0; q < 2; ++q) */ 18062fa5cd67SKarl Rupp if (n[0]) ++s; 1807af538404SDmitry Karpeev xstart += maxwidth; 1808af538404SDmitry Karpeev } /* for (i = 0; i < Mdomains; ++i) */ 18092bbb417fSDmitry Karpeev ystart += maxheight; 1810af538404SDmitry Karpeev } /* for (j = 0; j < Ndomains; ++j) */ 1811b1a0a8a3SJed Brown PetscFunctionReturn(0); 1812b1a0a8a3SJed Brown } 1813b1a0a8a3SJed Brown 1814b1a0a8a3SJed Brown /*@C 181506b43e7eSDmitry Karpeev PCGASMGetSubdomains - Gets the subdomains supported on this processor 181606b43e7eSDmitry Karpeev for the additive Schwarz preconditioner. 1817b1a0a8a3SJed Brown 1818b1a0a8a3SJed Brown Not Collective 1819b1a0a8a3SJed Brown 1820b1a0a8a3SJed Brown Input Parameter: 1821b1a0a8a3SJed Brown . pc - the preconditioner context 1822b1a0a8a3SJed Brown 1823b1a0a8a3SJed Brown Output Parameters: 1824b1a0a8a3SJed Brown + n - the number of subdomains for this processor (default value = 1) 18250298fd71SBarry Smith . iis - the index sets that define the inner subdomains (without overlap) supported on this processor (can be NULL) 18260298fd71SBarry Smith - ois - the index sets that define the outer subdomains (with overlap) supported on this processor (can be NULL) 1827b1a0a8a3SJed Brown 1828b1a0a8a3SJed Brown Notes: 18296a4f0f73SDmitry Karpeev The user is responsible for destroying the ISs and freeing the returned arrays. 1830b1a0a8a3SJed Brown The IS numbering is in the parallel, global numbering of the vector. 1831b1a0a8a3SJed Brown 1832b1a0a8a3SJed Brown Level: advanced 1833b1a0a8a3SJed Brown 1834db781477SPatrick Sanan .seealso: `PCGASMSetOverlap()`, `PCGASMGetSubKSP()`, `PCGASMCreateSubdomains2D()`, 1835db781477SPatrick Sanan `PCGASMSetSubdomains()`, `PCGASMGetSubmatrices()` 1836b1a0a8a3SJed Brown @*/ 18376a4f0f73SDmitry Karpeev PetscErrorCode PCGASMGetSubdomains(PC pc,PetscInt *n,IS *iis[],IS *ois[]) 1838b1a0a8a3SJed Brown { 1839f746d493SDmitry Karpeev PC_GASM *osm; 1840b1a0a8a3SJed Brown PetscBool match; 18416a4f0f73SDmitry Karpeev PetscInt i; 18425fd66863SKarl Rupp 1843b1a0a8a3SJed Brown PetscFunctionBegin; 1844b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 18459566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match)); 184628b400f6SJacob Faibussowitsch PetscCheck(match,PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Incorrect object type: expected %s, got %s instead", PCGASM, ((PetscObject)pc)->type_name); 1847f746d493SDmitry Karpeev osm = (PC_GASM*)pc->data; 1848ab3e923fSDmitry Karpeev if (n) *n = osm->n; 18496a4f0f73SDmitry Karpeev if (iis) { 18509566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(osm->n, iis)); 18516a4f0f73SDmitry Karpeev } 18526a4f0f73SDmitry Karpeev if (ois) { 18539566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(osm->n, ois)); 18546a4f0f73SDmitry Karpeev } 18556a4f0f73SDmitry Karpeev if (iis || ois) { 18566a4f0f73SDmitry Karpeev for (i = 0; i < osm->n; ++i) { 18576a4f0f73SDmitry Karpeev if (iis) (*iis)[i] = osm->iis[i]; 18586a4f0f73SDmitry Karpeev if (ois) (*ois)[i] = osm->ois[i]; 18596a4f0f73SDmitry Karpeev } 1860b1a0a8a3SJed Brown } 1861b1a0a8a3SJed Brown PetscFunctionReturn(0); 1862b1a0a8a3SJed Brown } 1863b1a0a8a3SJed Brown 1864b1a0a8a3SJed Brown /*@C 186506b43e7eSDmitry Karpeev PCGASMGetSubmatrices - Gets the local submatrices (for this processor 1866b1a0a8a3SJed Brown only) for the additive Schwarz preconditioner. 1867b1a0a8a3SJed Brown 1868b1a0a8a3SJed Brown Not Collective 1869b1a0a8a3SJed Brown 1870b1a0a8a3SJed Brown Input Parameter: 1871b1a0a8a3SJed Brown . pc - the preconditioner context 1872b1a0a8a3SJed Brown 1873b1a0a8a3SJed Brown Output Parameters: 1874b1a0a8a3SJed Brown + n - the number of matrices for this processor (default value = 1) 1875b1a0a8a3SJed Brown - mat - the matrices 1876b1a0a8a3SJed Brown 187795452b02SPatrick Sanan Notes: 187895452b02SPatrick Sanan matrices returned by this routine have the same communicators as the index sets (IS) 18798f3b4b4dSDmitry Karpeev used to define subdomains in PCGASMSetSubdomains() 1880b1a0a8a3SJed Brown Level: advanced 1881b1a0a8a3SJed Brown 1882db781477SPatrick Sanan .seealso: `PCGASMSetOverlap()`, `PCGASMGetSubKSP()`, 1883db781477SPatrick Sanan `PCGASMCreateSubdomains2D()`, `PCGASMSetSubdomains()`, `PCGASMGetSubdomains()` 1884b1a0a8a3SJed Brown @*/ 188506b43e7eSDmitry Karpeev PetscErrorCode PCGASMGetSubmatrices(PC pc,PetscInt *n,Mat *mat[]) 1886b1a0a8a3SJed Brown { 1887f746d493SDmitry Karpeev PC_GASM *osm; 1888b1a0a8a3SJed Brown PetscBool match; 1889b1a0a8a3SJed Brown 1890b1a0a8a3SJed Brown PetscFunctionBegin; 1891b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1892b1a0a8a3SJed Brown PetscValidIntPointer(n,2); 1893b1a0a8a3SJed Brown if (mat) PetscValidPointer(mat,3); 189428b400f6SJacob Faibussowitsch PetscCheck(pc->setupcalled,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call after KSPSetUp() or PCSetUp()."); 18959566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match)); 189628b400f6SJacob Faibussowitsch PetscCheck(match,PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Expected %s, got %s instead", PCGASM, ((PetscObject)pc)->type_name); 1897f746d493SDmitry Karpeev osm = (PC_GASM*)pc->data; 1898ab3e923fSDmitry Karpeev if (n) *n = osm->n; 1899b1a0a8a3SJed Brown if (mat) *mat = osm->pmat; 1900b1a0a8a3SJed Brown PetscFunctionReturn(0); 1901b1a0a8a3SJed Brown } 1902d709ea83SDmitry Karpeev 1903d709ea83SDmitry Karpeev /*@ 19048f3b4b4dSDmitry Karpeev PCGASMSetUseDMSubdomains - Indicates whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible. 1905d709ea83SDmitry Karpeev Logically Collective 1906d709ea83SDmitry Karpeev 1907d8d19677SJose E. Roman Input Parameters: 1908d709ea83SDmitry Karpeev + pc - the preconditioner 1909d709ea83SDmitry Karpeev - flg - boolean indicating whether to use subdomains defined by the DM 1910d709ea83SDmitry Karpeev 1911d709ea83SDmitry Karpeev Options Database Key: 19128f3b4b4dSDmitry Karpeev . -pc_gasm_dm_subdomains -pc_gasm_overlap -pc_gasm_total_subdomains 1913d709ea83SDmitry Karpeev 1914d709ea83SDmitry Karpeev Level: intermediate 1915d709ea83SDmitry Karpeev 1916d709ea83SDmitry Karpeev Notes: 19178f3b4b4dSDmitry Karpeev PCGASMSetSubdomains(), PCGASMSetTotalSubdomains() or PCGASMSetOverlap() take precedence over PCGASMSetUseDMSubdomains(), 19188f3b4b4dSDmitry Karpeev so setting PCGASMSetSubdomains() with nontrivial subdomain ISs or any of PCGASMSetTotalSubdomains() and PCGASMSetOverlap() 19198f3b4b4dSDmitry Karpeev automatically turns the latter off. 1920d709ea83SDmitry Karpeev 1921db781477SPatrick Sanan .seealso: `PCGASMGetUseDMSubdomains()`, `PCGASMSetSubdomains()`, `PCGASMSetOverlap()` 1922db781477SPatrick Sanan `PCGASMCreateSubdomains2D()` 1923d709ea83SDmitry Karpeev @*/ 19248f3b4b4dSDmitry Karpeev PetscErrorCode PCGASMSetUseDMSubdomains(PC pc,PetscBool flg) 1925d709ea83SDmitry Karpeev { 1926d709ea83SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 1927d709ea83SDmitry Karpeev PetscBool match; 1928d709ea83SDmitry Karpeev 1929d709ea83SDmitry Karpeev PetscFunctionBegin; 1930d709ea83SDmitry Karpeev PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1931d709ea83SDmitry Karpeev PetscValidLogicalCollectiveBool(pc,flg,2); 193228b400f6SJacob Faibussowitsch PetscCheck(!pc->setupcalled,((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Not for a setup PC."); 19339566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match)); 1934d709ea83SDmitry Karpeev if (match) { 19358f3b4b4dSDmitry Karpeev if (!osm->user_subdomains && osm->N == PETSC_DETERMINE && osm->overlap < 0) { 1936d709ea83SDmitry Karpeev osm->dm_subdomains = flg; 1937d709ea83SDmitry Karpeev } 19388f3b4b4dSDmitry Karpeev } 1939d709ea83SDmitry Karpeev PetscFunctionReturn(0); 1940d709ea83SDmitry Karpeev } 1941d709ea83SDmitry Karpeev 1942d709ea83SDmitry Karpeev /*@ 19438f3b4b4dSDmitry Karpeev PCGASMGetUseDMSubdomains - Returns flag indicating whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible. 1944d709ea83SDmitry Karpeev Not Collective 1945d709ea83SDmitry Karpeev 1946d709ea83SDmitry Karpeev Input Parameter: 1947d709ea83SDmitry Karpeev . pc - the preconditioner 1948d709ea83SDmitry Karpeev 1949d709ea83SDmitry Karpeev Output Parameter: 1950d709ea83SDmitry Karpeev . flg - boolean indicating whether to use subdomains defined by the DM 1951d709ea83SDmitry Karpeev 1952d709ea83SDmitry Karpeev Level: intermediate 1953d709ea83SDmitry Karpeev 1954db781477SPatrick Sanan .seealso: `PCGASMSetUseDMSubdomains()`, `PCGASMSetOverlap()` 1955db781477SPatrick Sanan `PCGASMCreateSubdomains2D()` 1956d709ea83SDmitry Karpeev @*/ 19578f3b4b4dSDmitry Karpeev PetscErrorCode PCGASMGetUseDMSubdomains(PC pc,PetscBool* flg) 1958d709ea83SDmitry Karpeev { 1959d709ea83SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 1960d709ea83SDmitry Karpeev PetscBool match; 1961d709ea83SDmitry Karpeev 1962d709ea83SDmitry Karpeev PetscFunctionBegin; 1963d709ea83SDmitry Karpeev PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1964534a8f05SLisandro Dalcin PetscValidBoolPointer(flg,2); 19659566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match)); 1966d709ea83SDmitry Karpeev if (match) { 1967d709ea83SDmitry Karpeev if (flg) *flg = osm->dm_subdomains; 1968d709ea83SDmitry Karpeev } 1969d709ea83SDmitry Karpeev PetscFunctionReturn(0); 1970d709ea83SDmitry Karpeev } 1971