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. */ 465f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc2(osm->n,numbering,osm->n,permutation)); 475f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectsListGetGlobalNumbering(PetscObjectComm((PetscObject)pc),osm->n,(PetscObject*)osm->iis,NULL,*numbering)); 488f3b4b4dSDmitry Karpeev for (i = 0; i < osm->n; ++i) (*permutation)[i] = i; 495f80ce2aSJacob Faibussowitsch CHKERRQ(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; 622c71b3e2SJacob Faibussowitsch PetscCheckFalse(i < 0 || i > osm->n,PetscObjectComm((PetscObject)viewer), PETSC_ERR_ARG_WRONG, "Invalid subdomain %D: must nonnegative and less than %D", i, osm->n); 634bde246dSDmitry Karpeev /* Inner subdomains. */ 645f80ce2aSJacob Faibussowitsch CHKERRQ(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 725f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(len, &cidx)); 735f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerStringOpen(PETSC_COMM_SELF, cidx, len, &sviewer)); 7436a9e3b9SBarry Smith #undef len 755f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetIndices(osm->iis[i], &idx)); 764bde246dSDmitry Karpeev for (j = 0; j < nidx; ++j) { 775f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerStringSPrintf(sviewer, "%D ", idx[j])); 784bde246dSDmitry Karpeev } 795f80ce2aSJacob Faibussowitsch CHKERRQ(ISRestoreIndices(osm->iis[i],&idx)); 805f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&sviewer)); 815f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer, "Inner subdomain:\n")); 825f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerFlush(viewer)); 835f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPushSynchronized(viewer)); 845f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIISynchronizedPrintf(viewer, "%s", cidx)); 855f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerFlush(viewer)); 865f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPopSynchronized(viewer)); 875f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer, "\n")); 885f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerFlush(viewer)); 895f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(cidx)); 904bde246dSDmitry Karpeev /* Outer subdomains. */ 915f80ce2aSJacob Faibussowitsch CHKERRQ(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 995f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(len, &cidx)); 1005f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerStringOpen(PETSC_COMM_SELF, cidx, len, &sviewer)); 10136a9e3b9SBarry Smith #undef len 1025f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetIndices(osm->ois[i], &idx)); 103af538404SDmitry Karpeev for (j = 0; j < nidx; ++j) { 1045f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerStringSPrintf(sviewer,"%D ", idx[j])); 105af538404SDmitry Karpeev } 1065f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerDestroy(&sviewer)); 1075f80ce2aSJacob Faibussowitsch CHKERRQ(ISRestoreIndices(osm->ois[i],&idx)); 1085f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer, "Outer subdomain:\n")); 1095f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerFlush(viewer)); 1105f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPushSynchronized(viewer)); 1115f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIISynchronizedPrintf(viewer, "%s", cidx)); 1125f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerFlush(viewer)); 1135f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPopSynchronized(viewer)); 1145f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer, "\n")); 1155f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerFlush(viewer)); 1165f80ce2aSJacob Faibussowitsch CHKERRQ(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; 1315f80ce2aSJacob Faibussowitsch CHKERRQ(PCGetOptionsPrefix(pc,&prefix)); 1325f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHasName(NULL,prefix,"-pc_gasm_print_subdomains",&found)); 1339140fbc5SPierre Jolivet if (!found) PetscFunctionReturn(0); 1345f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetString(NULL,prefix,"-pc_gasm_print_subdomains",fname,sizeof(fname),&found)); 1355f80ce2aSJacob Faibussowitsch if (!found) CHKERRQ(PetscStrcpy(fname,"stdout")); 1365f80ce2aSJacob Faibussowitsch CHKERRQ(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 */ 1415f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectName((PetscObject)viewer)); 1424bde246dSDmitry Karpeev l = 0; 1435f80ce2aSJacob Faibussowitsch CHKERRQ(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) { 1495f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerGetSubViewer(viewer,((PetscObject)osm->ois[d])->comm, &sviewer)); 1505f80ce2aSJacob Faibussowitsch CHKERRQ(PCGASMSubdomainView_Private(pc,d,sviewer)); 1515f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerRestoreSubViewer(viewer,((PetscObject)osm->ois[d])->comm, &sviewer)); 1524bde246dSDmitry Karpeev ++l; 153af538404SDmitry Karpeev } 1544bde246dSDmitry Karpeev } 1555f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Barrier(PetscObjectComm((PetscObject)pc))); 1564bde246dSDmitry Karpeev } 1575f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree2(numbering,permutation)); 1585f80ce2aSJacob Faibussowitsch CHKERRQ(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; 1775f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc), &size)); 1785f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc), &rank)); 17906b43e7eSDmitry Karpeev 18006b43e7eSDmitry Karpeev if (osm->overlap >= 0) { 1815f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSNPrintf(overlap,sizeof(overlap),"requested amount of overlap = %D",osm->overlap)); 18206b43e7eSDmitry Karpeev } 1838f3b4b4dSDmitry Karpeev if (osm->N != PETSC_DETERMINE) { 1845f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSNPrintf(gsubdomains, sizeof(gsubdomains), "total number of subdomains = %D",osm->N)); 18506b43e7eSDmitry Karpeev } 1868f3b4b4dSDmitry Karpeev if (osm->nmax != PETSC_DETERMINE) { 1875f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSNPrintf(msubdomains,sizeof(msubdomains),"max number of local subdomains = %D",osm->nmax)); 18806b43e7eSDmitry Karpeev } 18906b43e7eSDmitry Karpeev 1905f80ce2aSJacob Faibussowitsch CHKERRQ(PCGetOptionsPrefix(pc,&prefix)); 1915f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsGetBool(NULL,prefix,"-pc_gasm_view_subdomains",&view_subdomains,NULL)); 19206b43e7eSDmitry Karpeev 1935f80ce2aSJacob Faibussowitsch CHKERRQ(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 */ 2005f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectName((PetscObject)viewer)); 2015f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer," Restriction/interpolation type: %s\n",PCGASMTypes[osm->type])); 2025f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer," %s\n",overlap)); 2035f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer," %s\n",gsubdomains)); 2045f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer," %s\n",msubdomains)); 2055f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPushSynchronized(viewer)); 2065f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIISynchronizedPrintf(viewer," [%d|%d] number of locally-supported subdomains = %D\n",rank,size,osm->n)); 2075f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerFlush(viewer)); 2085f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPopSynchronized(viewer)); 2096a4f0f73SDmitry Karpeev /* Cannot take advantage of osm->same_subdomain_solvers without a global numbering of subdomains. */ 2105f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer," Subdomain solver info is as follows:\n")); 2115f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPushTab(viewer)); 2125f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(viewer," - - - - - - - - - - - - - - - - - -\n")); 21306b43e7eSDmitry Karpeev /* Make sure that everybody waits for the banner to be printed. */ 2145f80ce2aSJacob Faibussowitsch CHKERRMPI(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. */ 2165f80ce2aSJacob Faibussowitsch CHKERRQ(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) { 2235f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(((PetscObject)osm->ois[d])->comm, &ssize)); 2245f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(((PetscObject)osm->ois[d])->comm, &srank)); 2255f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerGetSubViewer(viewer,((PetscObject)osm->ois[d])->comm, &sviewer)); 2265f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetLocalSize(osm->ois[d],&bsz)); 2275f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIISynchronizedPrintf(sviewer," [%d|%d] (subcomm [%d|%d]) local subdomain number %D, local size = %D\n",rank,size,srank,ssize,d,bsz)); 2285f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerFlush(sviewer)); 22906b43e7eSDmitry Karpeev if (view_subdomains) { 2305f80ce2aSJacob Faibussowitsch CHKERRQ(PCGASMSubdomainView_Private(pc,d,sviewer)); 23106b43e7eSDmitry Karpeev } 2326a4f0f73SDmitry Karpeev if (!pc->setupcalled) { 2335f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(sviewer, " Solver not set up yet: PCSetUp() not yet called\n")); 2342fa5cd67SKarl Rupp } else { 2355f80ce2aSJacob Faibussowitsch CHKERRQ(KSPView(osm->ksp[d],sviewer)); 2366a4f0f73SDmitry Karpeev } 2375f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPrintf(sviewer," - - - - - - - - - - - - - - - - - -\n")); 2385f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerFlush(sviewer)); 2395f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerRestoreSubViewer(viewer,((PetscObject)osm->ois[d])->comm, &sviewer)); 24006b43e7eSDmitry Karpeev ++l; 241b1a0a8a3SJed Brown } 24206b43e7eSDmitry Karpeev } 2435f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Barrier(PetscObjectComm((PetscObject)pc))); 244b1a0a8a3SJed Brown } 2455f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree2(numbering,permutation)); 2465f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerASCIIPopTab(viewer)); 2475f80ce2aSJacob Faibussowitsch CHKERRQ(PetscViewerFlush(viewer)); 2489530cbd7SBarry Smith /* this line is needed to match the extra PetscViewerASCIIPushSynchronized() in PetscViewerGetSubViewer() */ 2495f80ce2aSJacob Faibussowitsch CHKERRQ(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; 2675f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectGetComm((PetscObject)pc,&comm)); 2685f80ce2aSJacob Faibussowitsch CHKERRMPI(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); 2752c71b3e2SJacob Faibussowitsch PetscCheckFalse(size%osm->N != 0,PETSC_COMM_WORLD,PETSC_ERR_ARG_INCOMP,"have to specify the total number of subdomains %D to be a factor of the number of processors %d ",osm->N,size); 276ea91fabdSFande Kong nlocalsubdomains = size/osm->N; 277ea91fabdSFande Kong osm->n = 1; 2785f80ce2aSJacob Faibussowitsch CHKERRQ(MatPartitioningCreate(comm,&part)); 2795f80ce2aSJacob Faibussowitsch CHKERRQ(MatPartitioningSetAdjacency(part,pc->pmat)); 2805f80ce2aSJacob Faibussowitsch CHKERRQ(MatPartitioningSetType(part,MATPARTITIONINGHIERARCH)); 2815f80ce2aSJacob Faibussowitsch CHKERRQ(MatPartitioningHierarchicalSetNcoarseparts(part,osm->N)); 2825f80ce2aSJacob Faibussowitsch CHKERRQ(MatPartitioningHierarchicalSetNfineparts(part,nlocalsubdomains)); 2835f80ce2aSJacob Faibussowitsch CHKERRQ(MatPartitioningSetFromOptions(part)); 284ea91fabdSFande Kong /* get new processor owner number of each vertex */ 2855f80ce2aSJacob Faibussowitsch CHKERRQ(MatPartitioningApply(part,&partitioning)); 2865f80ce2aSJacob Faibussowitsch CHKERRQ(ISBuildTwoSided(partitioning,NULL,&fromrows)); 2875f80ce2aSJacob Faibussowitsch CHKERRQ(ISPartitioningToNumbering(partitioning,&isn)); 2885f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&isn)); 2895f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetLocalSize(fromrows,&fromrows_localsize)); 2905f80ce2aSJacob Faibussowitsch CHKERRQ(MatPartitioningDestroy(&part)); 2915f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(pc->pmat,&outervec,NULL)); 2925f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateMPI(comm,fromrows_localsize,PETSC_DETERMINE,&(osm->pcx))); 2935f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(osm->pcx,&(osm->pcy))); 2945f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterCreate(osm->pcx,NULL,outervec,fromrows,&(osm->pctoouter))); 2955f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSubMatrix(pc->pmat,fromrows,fromrows,MAT_INITIAL_MATRIX,&(osm->permutationP))); 2965f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectReference((PetscObject)fromrows)); 297ea91fabdSFande Kong osm->permutationIS = fromrows; 298ea91fabdSFande Kong osm->pcmat = pc->pmat; 2995f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectReference((PetscObject)osm->permutationP)); 300ea91fabdSFande Kong pc->pmat = osm->permutationP; 3015f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&outervec)); 3025f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&fromrows)); 3035f80ce2aSJacob Faibussowitsch CHKERRQ(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; 3325f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size)); 3335f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank)); 334b1a0a8a3SJed Brown if (!pc->setupcalled) { 335ea91fabdSFande Kong /* use a hierarchical partitioning */ 336ea91fabdSFande Kong if (osm->hierarchicalpartitioning) { 3375f80ce2aSJacob Faibussowitsch CHKERRQ(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. */ 3425f80ce2aSJacob Faibussowitsch CHKERRQ(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; 3475f80ce2aSJacob Faibussowitsch CHKERRQ(DMCreateDomainDecomposition(pc->dm, &num_subdomains, &subdomain_names, &inner_subdomain_is, &outer_subdomain_is, &subdomain_dm)); 348a35b7b57SDmitry Karpeev if (num_subdomains) { 3495f80ce2aSJacob Faibussowitsch CHKERRQ(PCGASMSetSubdomains(pc, num_subdomains, inner_subdomain_is, outer_subdomain_is)); 35069ca1f37SDmitry Karpeev } 351a35b7b57SDmitry Karpeev for (d = 0; d < num_subdomains; ++d) { 3525f80ce2aSJacob Faibussowitsch if (inner_subdomain_is) CHKERRQ(ISDestroy(&inner_subdomain_is[d])); 3535f80ce2aSJacob Faibussowitsch if (outer_subdomain_is) CHKERRQ(ISDestroy(&outer_subdomain_is[d])); 354fdc48646SDmitry Karpeev } 3555f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(inner_subdomain_is)); 3565f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(outer_subdomain_is)); 3578f3b4b4dSDmitry Karpeev } else { 3588f3b4b4dSDmitry Karpeev /* still no subdomains; use one per processor */ 35944fe18e5SDmitry Karpeev osm->nmax = osm->n = 1; 3605f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size)); 361f746d493SDmitry Karpeev osm->N = size; 3625f80ce2aSJacob Faibussowitsch CHKERRQ(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 */ 3715f80ce2aSJacob Faibussowitsch CHKERRQ(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 */ 3785f80ce2aSJacob Faibussowitsch CHKERRQ(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 */ 3815f80ce2aSJacob Faibussowitsch CHKERRQ(ISDuplicate(osm->iis[i],(osm->ois)+i)); 3825f80ce2aSJacob Faibussowitsch CHKERRQ(ISCopy(osm->iis[i],osm->ois[i])); 3838f3b4b4dSDmitry Karpeev } else { 3845f80ce2aSJacob Faibussowitsch CHKERRQ(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 */ 3905f80ce2aSJacob Faibussowitsch CHKERRQ(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; 3995f80ce2aSJacob Faibussowitsch CHKERRMPI(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. */ 4045f80ce2aSJacob Faibussowitsch CHKERRQ(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++) { 4095f80ce2aSJacob Faibussowitsch CHKERRQ(ISSort(osm->ois[i])); 4105f80ce2aSJacob Faibussowitsch CHKERRQ(ISSort(osm->iis[i])); 411b1a0a8a3SJed Brown } 412b1a0a8a3SJed Brown } 4135f80ce2aSJacob Faibussowitsch CHKERRQ(PCGetOptionsPrefix(pc,&prefix)); 4145f80ce2aSJacob Faibussowitsch CHKERRQ(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++) { 4225f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetLocalSize(osm->ois[i],&oni)); 4236a4f0f73SDmitry Karpeev on += oni; 424f746d493SDmitry Karpeev } 4255f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(on, &oidx)); 4266a4f0f73SDmitry Karpeev on = 0; 427930d09c1SFande Kong /* Merge local indices together */ 428f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 4295f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetLocalSize(osm->ois[i],&oni)); 4305f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetIndices(osm->ois[i],&oidxi)); 4315f80ce2aSJacob Faibussowitsch CHKERRQ(PetscArraycpy(oidx+on,oidxi,oni)); 4325f80ce2aSJacob Faibussowitsch CHKERRQ(ISRestoreIndices(osm->ois[i],&oidxi)); 4336a4f0f73SDmitry Karpeev on += oni; 434f746d493SDmitry Karpeev } 4355f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateGeneral(((PetscObject)(pc))->comm,on,oidx,PETSC_OWN_POINTER,&gois)); 436ea91fabdSFande Kong nTotalInnerIndices = 0; 437ea91fabdSFande Kong for (i=0; i<osm->n; i++) { 4385f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetLocalSize(osm->iis[i],&nInnerIndices)); 439ea91fabdSFande Kong nTotalInnerIndices += nInnerIndices; 440ea91fabdSFande Kong } 4415f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateMPI(((PetscObject)(pc))->comm,nTotalInnerIndices,PETSC_DETERMINE,&x)); 4425f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(x,&y)); 443ea91fabdSFande Kong 4445f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateMPI(PetscObjectComm((PetscObject)pc),on,PETSC_DECIDE,&osm->gx)); 4455f80ce2aSJacob Faibussowitsch CHKERRQ(VecDuplicate(osm->gx,&osm->gy)); 4465f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetOwnershipRange(osm->gx, &gostart, NULL)); 4475f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateStride(PetscObjectComm((PetscObject)pc),on,gostart,1, &goid)); 448930d09c1SFande Kong /* gois might indices not on local */ 4495f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterCreate(x,gois,osm->gx,goid, &(osm->gorestriction))); 4505f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(osm->n,&numbering)); 4515f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectsListGetGlobalNumbering(PetscObjectComm((PetscObject)pc),osm->n,(PetscObject*)osm->ois,NULL,numbering)); 4525f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&x)); 4535f80ce2aSJacob Faibussowitsch CHKERRQ(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++) { 4685f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetLocalSize(osm->ois[i],&oni)); 4696a4f0f73SDmitry Karpeev on += oni; 470f746d493SDmitry Karpeev } 4715f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(on, &iidx)); 4725f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(on, &ioidx)); 4735f80ce2aSJacob Faibussowitsch CHKERRQ(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++) { 4775f80ce2aSJacob Faibussowitsch CHKERRQ(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 } 4835f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(y,&array)); 4845f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(osm->gorestriction,y,osm->gy,INSERT_VALUES,SCATTER_FORWARD)); 4855f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(osm->gorestriction,y,osm->gy,INSERT_VALUES,SCATTER_FORWARD)); 4865f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetOwnershipRange(osm->gy,&gostart, NULL)); 4875f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(osm->gy,&array)); 488f771a274SFande Kong on = 0; 489f771a274SFande Kong in = 0; 490f771a274SFande Kong for (i=0; i<osm->n; i++) { 4915f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetLocalSize(osm->ois[i],&oni)); 4925f80ce2aSJacob Faibussowitsch CHKERRQ(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 } 4995f80ce2aSJacob Faibussowitsch CHKERRQ(ISRestoreIndices(osm->ois[i], &indices)); 500f771a274SFande Kong on += oni; 501f771a274SFande Kong } 5025f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(osm->gy,&array)); 5035f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateGeneral(PetscObjectComm((PetscObject)pc),in,iidx,PETSC_OWN_POINTER,&giis)); 5045f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateGeneral(PetscObjectComm((PetscObject)pc),in,ioidx,PETSC_OWN_POINTER,&giois)); 5055f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterCreate(y,giis,osm->gy,giois,&osm->girestriction)); 5065f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&y)); 5075f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&giis)); 5085f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&giois)); 509b1a0a8a3SJed Brown } 5105f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&goid)); 5115f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(numbering)); 5122fa5cd67SKarl Rupp 513f746d493SDmitry Karpeev /* Create the subdomain work vectors. */ 5145f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(osm->n,&osm->x)); 5155f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(osm->n,&osm->y)); 5165f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(osm->gx, &gxarray)); 5175f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(osm->gy, &gyarray)); 5186a4f0f73SDmitry Karpeev for (i=0, on=0; i<osm->n; ++i, on += oni) { 5196a4f0f73SDmitry Karpeev PetscInt oNi; 5205f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetLocalSize(osm->ois[i],&oni)); 521930d09c1SFande Kong /* on a sub communicator */ 5225f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetSize(osm->ois[i],&oNi)); 5235f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateMPIWithArray(((PetscObject)(osm->ois[i]))->comm,1,oni,oNi,gxarray+on,&osm->x[i])); 5245f80ce2aSJacob Faibussowitsch CHKERRQ(VecCreateMPIWithArray(((PetscObject)(osm->ois[i]))->comm,1,oni,oNi,gyarray+on,&osm->y[i])); 525b1a0a8a3SJed Brown } 5265f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(osm->gx, &gxarray)); 5275f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(osm->gy, &gyarray)); 5288f3b4b4dSDmitry Karpeev /* Create the subdomain solvers */ 5295f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(osm->n,&osm->ksp)); 5308f3b4b4dSDmitry Karpeev for (i=0; i<osm->n; i++) { 5318f3b4b4dSDmitry Karpeev char subprefix[PETSC_MAX_PATH_LEN+1]; 5325f80ce2aSJacob Faibussowitsch CHKERRQ(KSPCreate(((PetscObject)(osm->ois[i]))->comm,&ksp)); 5335f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSetErrorIfNotConverged(ksp,pc->erroriffailure)); 5345f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogObjectParent((PetscObject)pc,(PetscObject)ksp)); 5355f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectIncrementTabLevel((PetscObject)ksp,(PetscObject)pc,1)); 5365f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSetType(ksp,KSPPREONLY)); 5375f80ce2aSJacob Faibussowitsch CHKERRQ(KSPGetPC(ksp,&subpc)); /* Why do we need this here? */ 5388f3b4b4dSDmitry Karpeev if (subdomain_dm) { 5395f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSetDM(ksp,subdomain_dm[i])); 5405f80ce2aSJacob Faibussowitsch CHKERRQ(DMDestroy(subdomain_dm+i)); 5418f3b4b4dSDmitry Karpeev } 5425f80ce2aSJacob Faibussowitsch CHKERRQ(PCGetOptionsPrefix(pc,&prefix)); 5435f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSetOptionsPrefix(ksp,prefix)); 5448f3b4b4dSDmitry Karpeev if (subdomain_names && subdomain_names[i]) { 5455f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSNPrintf(subprefix,PETSC_MAX_PATH_LEN,"sub_%s_",subdomain_names[i])); 5465f80ce2aSJacob Faibussowitsch CHKERRQ(KSPAppendOptionsPrefix(ksp,subprefix)); 5475f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(subdomain_names[i])); 5488f3b4b4dSDmitry Karpeev } 5495f80ce2aSJacob Faibussowitsch CHKERRQ(KSPAppendOptionsPrefix(ksp,"sub_")); 550b1a0a8a3SJed Brown osm->ksp[i] = ksp; 551b1a0a8a3SJed Brown } 5525f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(subdomain_dm)); 5535f80ce2aSJacob Faibussowitsch CHKERRQ(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) { 5605f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroyMatrices(osm->n,&osm->pmat)); 561b1a0a8a3SJed Brown scall = MAT_INITIAL_MATRIX; 562b1a0a8a3SJed Brown } 563ea91fabdSFande Kong if (osm->permutationIS) { 5645f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSubMatrix(pc->pmat,osm->permutationIS,osm->permutationIS,scall,&osm->permutationP)); 5655f80ce2aSJacob Faibussowitsch CHKERRQ(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) { 5755f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSubMatricesMPI(pc->pmat,osm->n,osm->ois,osm->ois,scall,&osm->pmat)); 5762fa5cd67SKarl Rupp } else { 5775f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateSubMatrices(pc->pmat,osm->n,osm->ois,osm->ois,scall,&osm->pmat)); 57881a5c4aaSDmitry Karpeev } 579b1a0a8a3SJed Brown if (scall == MAT_INITIAL_MATRIX) { 5805f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectGetOptionsPrefix((PetscObject)pc->pmat,&pprefix)); 581f746d493SDmitry Karpeev for (i=0; i<osm->n; i++) { 5825f80ce2aSJacob Faibussowitsch CHKERRQ(PetscLogObjectParent((PetscObject)pc,(PetscObject)osm->pmat[i])); 5835f80ce2aSJacob Faibussowitsch CHKERRQ(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) */ 5895f80ce2aSJacob Faibussowitsch CHKERRQ(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++) { 5955f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSetOperators(osm->ksp[i],osm->pmat[i],osm->pmat[i])); 5965f80ce2aSJacob Faibussowitsch CHKERRQ(KSPGetOptionsPrefix(osm->ksp[i],&prefix)); 5975f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOptionsPrefix(osm->pmat[i],prefix)); 598b1a0a8a3SJed Brown if (!pc->setupcalled) { 5995f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSetFromOptions(osm->ksp[i])); 600b1a0a8a3SJed Brown } 601b1a0a8a3SJed Brown } 602ea91fabdSFande Kong if (osm->pcmat) { 6035f80ce2aSJacob Faibussowitsch CHKERRQ(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++) { 6175f80ce2aSJacob Faibussowitsch CHKERRQ(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) { 6315f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(osm->pctoouter,xin,osm->pcx,INSERT_VALUES,SCATTER_REVERSE)); 6325f80ce2aSJacob Faibussowitsch CHKERRQ(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 */ 6455f80ce2aSJacob Faibussowitsch CHKERRQ(VecZeroEntries(osm->gx)); 6465f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(osm->girestriction,x,osm->gx,INSERT_VALUES,forward)); 6472fa5cd67SKarl Rupp } else { 6485f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward)); 649b1a0a8a3SJed Brown } 6505f80ce2aSJacob Faibussowitsch CHKERRQ(VecZeroEntries(osm->gy)); 6516a4f0f73SDmitry Karpeev if (!(osm->type & PC_GASM_RESTRICT)) { 6525f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(osm->girestriction,x,osm->gx,INSERT_VALUES,forward)); 6532fa5cd67SKarl Rupp } else { 6545f80ce2aSJacob Faibussowitsch CHKERRQ(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) { 6585f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSolve(osm->ksp[i],osm->x[i],osm->y[i])); 6595f80ce2aSJacob Faibussowitsch CHKERRQ(KSPCheckSolve(osm->ksp[i],pc,osm->y[i])); 660b1a0a8a3SJed Brown } 66148e38a8aSPierre Jolivet /* do we need to zero y? */ 6625f80ce2aSJacob Faibussowitsch CHKERRQ(VecZeroEntries(y)); 6636a4f0f73SDmitry Karpeev if (!(osm->type & PC_GASM_INTERPOLATE)) { 6645f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(osm->girestriction,osm->gy,y,ADD_VALUES,reverse)); 6655f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(osm->girestriction,osm->gy,y,ADD_VALUES,reverse)); 6662fa5cd67SKarl Rupp } else { 6675f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse)); 6685f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse)); 6696a4f0f73SDmitry Karpeev } 670ea91fabdSFande Kong if (osm->pctoouter) { 6715f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(osm->pctoouter,y,yout,INSERT_VALUES,SCATTER_FORWARD)); 6725f80ce2aSJacob Faibussowitsch CHKERRQ(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; 6862c71b3e2SJacob Faibussowitsch PetscCheckFalse(osm->n != 1,PetscObjectComm((PetscObject)pc),PETSC_ERR_SUP,"Not yet implemented"); 6875f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetSize(Xin,NULL,&N)); 68848e38a8aSPierre Jolivet if (osm->pctoouter) { 6895f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetLocalSize(osm->pcx,&m)); 6905f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetSize(osm->pcx,&M)); 6915f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateDense(PetscObjectComm((PetscObject)osm->ois[0]),m,PETSC_DECIDE,M,N,NULL,&O)); 69248e38a8aSPierre Jolivet for (i = 0; i < N; ++i) { 6935f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseGetColumnVecRead(Xin,i,&x)); 6945f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseGetColumnVecWrite(O,i,&y)); 6955f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(osm->pctoouter,x,y,INSERT_VALUES,SCATTER_REVERSE)); 6965f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(osm->pctoouter,x,y,INSERT_VALUES,SCATTER_REVERSE)); 6975f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseRestoreColumnVecWrite(O,i,&y)); 6985f80ce2aSJacob Faibussowitsch CHKERRQ(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 */ 7095f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetLocalSize(osm->x[0],&m)); 7105f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetSize(osm->x[0],&M)); 7115f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateDense(PetscObjectComm((PetscObject)osm->ois[0]),m,PETSC_DECIDE,M,N,NULL,&Z)); 71248e38a8aSPierre Jolivet for (i = 0; i < N; ++i) { 7135f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseGetColumnVecRead(X,i,&x)); 7145f80ce2aSJacob Faibussowitsch CHKERRQ(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 */ 7175f80ce2aSJacob Faibussowitsch CHKERRQ(VecZeroEntries(y)); 7185f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(osm->girestriction,x,y,INSERT_VALUES,forward)); 7195f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(osm->girestriction,x,y,INSERT_VALUES,forward)); 72048e38a8aSPierre Jolivet } else { 7215f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(osm->gorestriction,x,y,INSERT_VALUES,forward)); 7225f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(osm->gorestriction,x,y,INSERT_VALUES,forward)); 72348e38a8aSPierre Jolivet } 7245f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseRestoreColumnVecWrite(Z,i,&y)); 7255f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseRestoreColumnVecRead(X,i,&x)); 72648e38a8aSPierre Jolivet } 7275f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateDense(PetscObjectComm((PetscObject)osm->ois[0]),m,PETSC_DECIDE,M,N,NULL,&W)); 7285f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetOption(Z,MAT_NO_OFF_PROC_ENTRIES,PETSC_TRUE)); 7295f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(Z,MAT_FINAL_ASSEMBLY)); 7305f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(Z,MAT_FINAL_ASSEMBLY)); 73148e38a8aSPierre Jolivet /* do the subdomain solve */ 7325f80ce2aSJacob Faibussowitsch CHKERRQ(KSPMatSolve(osm->ksp[0],Z,W)); 7335f80ce2aSJacob Faibussowitsch CHKERRQ(KSPCheckSolve(osm->ksp[0],pc,NULL)); 7345f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&Z)); 73548e38a8aSPierre Jolivet /* do we need to zero y? */ 7365f80ce2aSJacob Faibussowitsch CHKERRQ(MatZeroEntries(Y)); 73748e38a8aSPierre Jolivet for (i = 0; i < N; ++i) { 7385f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseGetColumnVecWrite(Y,i,&y)); 7395f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseGetColumnVecRead(W,i,&x)); 74048e38a8aSPierre Jolivet if (!(osm->type & PC_GASM_INTERPOLATE)) { 7415f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(osm->girestriction,x,y,ADD_VALUES,reverse)); 7425f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(osm->girestriction,x,y,ADD_VALUES,reverse)); 74348e38a8aSPierre Jolivet } else { 7445f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(osm->gorestriction,x,y,ADD_VALUES,reverse)); 7455f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(osm->gorestriction,x,y,ADD_VALUES,reverse)); 74648e38a8aSPierre Jolivet } 7475f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseRestoreColumnVecRead(W,i,&x)); 74848e38a8aSPierre Jolivet if (osm->pctoouter) { 7495f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseGetColumnVecWrite(Yout,i,&x)); 7505f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(osm->pctoouter,y,x,INSERT_VALUES,SCATTER_FORWARD)); 7515f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(osm->pctoouter,y,x,INSERT_VALUES,SCATTER_FORWARD)); 7525f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseRestoreColumnVecRead(Yout,i,&x)); 75348e38a8aSPierre Jolivet } 7545f80ce2aSJacob Faibussowitsch CHKERRQ(MatDenseRestoreColumnVecWrite(Y,i,&y)); 75548e38a8aSPierre Jolivet } 7565f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&W)); 7575f80ce2aSJacob Faibussowitsch CHKERRQ(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) { 7705f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(osm->pctoouter,xin,osm->pcx,INSERT_VALUES,SCATTER_REVERSE)); 7715f80ce2aSJacob Faibussowitsch CHKERRQ(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 */ 7875f80ce2aSJacob Faibussowitsch CHKERRQ(VecZeroEntries(osm->gx)); 7885f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(osm->girestriction,x,osm->gx,INSERT_VALUES,forward)); 7892fa5cd67SKarl Rupp } else { 7905f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(osm->gorestriction,x,osm->gx,INSERT_VALUES,forward)); 791b1a0a8a3SJed Brown } 7925f80ce2aSJacob Faibussowitsch CHKERRQ(VecZeroEntries(osm->gy)); 7936a4f0f73SDmitry Karpeev if (!(osm->type & PC_GASM_INTERPOLATE)) { 7945f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(osm->girestriction,x,osm->gx,INSERT_VALUES,forward)); 7952fa5cd67SKarl Rupp } else { 7965f80ce2aSJacob Faibussowitsch CHKERRQ(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. */ 8005f80ce2aSJacob Faibussowitsch CHKERRQ(KSPSolveTranspose(osm->ksp[i],osm->x[i],osm->y[i])); 8015f80ce2aSJacob Faibussowitsch CHKERRQ(KSPCheckSolve(osm->ksp[i],pc,osm->y[i])); 802b1a0a8a3SJed Brown } 8035f80ce2aSJacob Faibussowitsch CHKERRQ(VecZeroEntries(y)); 8046a4f0f73SDmitry Karpeev if (!(osm->type & PC_GASM_RESTRICT)) { 8055f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(osm->girestriction,osm->gy,y,ADD_VALUES,reverse)); 8065f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(osm->girestriction,osm->gy,y,ADD_VALUES,reverse)); 8072fa5cd67SKarl Rupp } else { 8085f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse)); 8095f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterEnd(osm->gorestriction,osm->gy,y,ADD_VALUES,reverse)); 8106a4f0f73SDmitry Karpeev } 811ea91fabdSFande Kong if (osm->pctoouter) { 8125f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterBegin(osm->pctoouter,y,yout,INSERT_VALUES,SCATTER_FORWARD)); 8135f80ce2aSJacob Faibussowitsch CHKERRQ(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++) { 8265f80ce2aSJacob Faibussowitsch CHKERRQ(KSPReset(osm->ksp[i])); 827b1a0a8a3SJed Brown } 828b1a0a8a3SJed Brown } 829b1a0a8a3SJed Brown if (osm->pmat) { 830f746d493SDmitry Karpeev if (osm->n > 0) { 831df750dc8SHong Zhang PetscMPIInt size; 8325f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size)); 833df750dc8SHong Zhang if (size > 1) { 8347dae84e0SHong Zhang /* osm->pmat is created by MatCreateSubMatricesMPI(), cannot use MatDestroySubMatrices() */ 8355f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroyMatrices(osm->n,&osm->pmat)); 836df750dc8SHong Zhang } else { 8375f80ce2aSJacob Faibussowitsch CHKERRQ(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++) { 8435f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&osm->x[i])); 8445f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&osm->y[i])); 845b1a0a8a3SJed Brown } 846d34fcf5fSBarry Smith } 8475f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&osm->gx)); 8485f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&osm->gy)); 849ab3e923fSDmitry Karpeev 8505f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterDestroy(&osm->gorestriction)); 8515f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterDestroy(&osm->girestriction)); 8528f3b4b4dSDmitry Karpeev if (!osm->user_subdomains) { 8535f80ce2aSJacob Faibussowitsch CHKERRQ(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) { 8585f80ce2aSJacob Faibussowitsch CHKERRQ(VecScatterDestroy(&(osm->pctoouter))); 859ea91fabdSFande Kong } 860ea91fabdSFande Kong if (osm->permutationIS) { 8615f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&(osm->permutationIS))); 862ea91fabdSFande Kong } 863ea91fabdSFande Kong if (osm->pcx) { 8645f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&(osm->pcx))); 865ea91fabdSFande Kong } 866ea91fabdSFande Kong if (osm->pcy) { 8675f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&(osm->pcy))); 868ea91fabdSFande Kong } 869ea91fabdSFande Kong if (osm->permutationP) { 8705f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&(osm->permutationP))); 871ea91fabdSFande Kong } 872ea91fabdSFande Kong if (osm->pcmat) { 8735f80ce2aSJacob Faibussowitsch CHKERRQ(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; 8845f80ce2aSJacob Faibussowitsch CHKERRQ(PCReset_GASM(pc)); 885135757f6SDmitry Karpeev /* PCReset will not destroy subdomains, if user_subdomains is true. */ 8865f80ce2aSJacob Faibussowitsch CHKERRQ(PCGASMDestroySubdomains(osm->n,&osm->ois,&osm->iis)); 887a06653b4SBarry Smith if (osm->ksp) { 888a06653b4SBarry Smith for (i=0; i<osm->n; i++) { 8895f80ce2aSJacob Faibussowitsch CHKERRQ(KSPDestroy(&osm->ksp[i])); 890a06653b4SBarry Smith } 8915f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(osm->ksp)); 892a06653b4SBarry Smith } 8935f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(osm->x)); 8945f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(osm->y)); 8955f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(pc->data)); 896b1a0a8a3SJed Brown PetscFunctionReturn(0); 897b1a0a8a3SJed Brown } 898b1a0a8a3SJed Brown 8994416b707SBarry Smith static PetscErrorCode PCSetFromOptions_GASM(PetscOptionItems *PetscOptionsObject,PC pc) 900a6dfd86eSKarl Rupp { 901f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 902b1a0a8a3SJed Brown PetscInt blocks,ovl; 90357501b6eSBarry Smith PetscBool flg; 904f746d493SDmitry Karpeev PCGASMType gasmtype; 905b1a0a8a3SJed Brown 906b1a0a8a3SJed Brown PetscFunctionBegin; 9075f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsHead(PetscOptionsObject,"Generalized additive Schwarz options")); 9085f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-pc_gasm_use_dm_subdomains","If subdomains aren't set, use DMCreateDomainDecomposition() to define subdomains.","PCGASMSetUseDMSubdomains",osm->dm_subdomains,&osm->dm_subdomains,&flg)); 9095f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsInt("-pc_gasm_total_subdomains","Total number of subdomains across communicator","PCGASMSetTotalSubdomains",osm->N,&blocks,&flg)); 91065db9045SDmitry Karpeev if (flg) { 9115f80ce2aSJacob Faibussowitsch CHKERRQ(PCGASMSetTotalSubdomains(pc,blocks)); 91265db9045SDmitry Karpeev } 9135f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsInt("-pc_gasm_overlap","Number of overlapping degrees of freedom","PCGASMSetOverlap",osm->overlap,&ovl,&flg)); 91465db9045SDmitry Karpeev if (flg) { 9155f80ce2aSJacob Faibussowitsch CHKERRQ(PCGASMSetOverlap(pc,ovl)); 916d709ea83SDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 91765db9045SDmitry Karpeev } 918b1a0a8a3SJed Brown flg = PETSC_FALSE; 9195f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsEnum("-pc_gasm_type","Type of restriction/extension","PCGASMSetType",PCGASMTypes,(PetscEnum)osm->type,(PetscEnum*)&gasmtype,&flg)); 9205f80ce2aSJacob Faibussowitsch if (flg) CHKERRQ(PCGASMSetType(pc,gasmtype)); 9215f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-pc_gasm_use_hierachical_partitioning","use hierarchical partitioning",NULL,osm->hierarchicalpartitioning,&osm->hierarchicalpartitioning,&flg)); 9225f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsTail()); 923b1a0a8a3SJed Brown PetscFunctionReturn(0); 924b1a0a8a3SJed Brown } 925b1a0a8a3SJed Brown 926b1a0a8a3SJed Brown /*------------------------------------------------------------------------------------*/ 927b1a0a8a3SJed Brown 9288f3b4b4dSDmitry Karpeev /*@ 9298f3b4b4dSDmitry Karpeev PCGASMSetTotalSubdomains - sets the total number of subdomains to use across the 9308f3b4b4dSDmitry Karpeev communicator. 9318f3b4b4dSDmitry Karpeev Logically collective on pc 9328f3b4b4dSDmitry Karpeev 9338f3b4b4dSDmitry Karpeev Input Parameters: 9348f3b4b4dSDmitry Karpeev + pc - the preconditioner 9358f3b4b4dSDmitry Karpeev - N - total number of subdomains 9368f3b4b4dSDmitry Karpeev 9378f3b4b4dSDmitry Karpeev Level: beginner 9388f3b4b4dSDmitry Karpeev 9398f3b4b4dSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMSetOverlap() 9408f3b4b4dSDmitry Karpeev PCGASMCreateSubdomains2D() 9418f3b4b4dSDmitry Karpeev @*/ 9428f3b4b4dSDmitry Karpeev PetscErrorCode PCGASMSetTotalSubdomains(PC pc,PetscInt N) 9438f3b4b4dSDmitry Karpeev { 9448f3b4b4dSDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 9458f3b4b4dSDmitry Karpeev PetscMPIInt size,rank; 9468f3b4b4dSDmitry Karpeev 9478f3b4b4dSDmitry Karpeev PetscFunctionBegin; 9482c71b3e2SJacob Faibussowitsch PetscCheckFalse(N < 1,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Total number of subdomains must be 1 or more, got N = %D",N); 949*28b400f6SJacob Faibussowitsch PetscCheck(!pc->setupcalled,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetTotalSubdomains() should be called before calling PCSetUp()."); 9508f3b4b4dSDmitry Karpeev 9515f80ce2aSJacob Faibussowitsch CHKERRQ(PCGASMDestroySubdomains(osm->n,&osm->iis,&osm->ois)); 9528f3b4b4dSDmitry Karpeev osm->ois = osm->iis = NULL; 9538f3b4b4dSDmitry Karpeev 9545f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size)); 9555f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)pc),&rank)); 9568f3b4b4dSDmitry Karpeev osm->N = N; 9578f3b4b4dSDmitry Karpeev osm->n = PETSC_DETERMINE; 9588f3b4b4dSDmitry Karpeev osm->nmax = PETSC_DETERMINE; 9598f3b4b4dSDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 9608f3b4b4dSDmitry Karpeev PetscFunctionReturn(0); 9618f3b4b4dSDmitry Karpeev } 9628f3b4b4dSDmitry Karpeev 963b541e6a4SDmitry Karpeev static PetscErrorCode PCGASMSetSubdomains_GASM(PC pc,PetscInt n,IS iis[],IS ois[]) 964b1a0a8a3SJed Brown { 965f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 966b1a0a8a3SJed Brown PetscInt i; 967b1a0a8a3SJed Brown 968b1a0a8a3SJed Brown PetscFunctionBegin; 9692c71b3e2SJacob Faibussowitsch PetscCheckFalse(n < 1,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Each process must have 1 or more subdomains, got n = %D",n); 970*28b400f6SJacob Faibussowitsch PetscCheck(!pc->setupcalled,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetSubdomains() should be called before calling PCSetUp()."); 971b1a0a8a3SJed Brown 9725f80ce2aSJacob Faibussowitsch CHKERRQ(PCGASMDestroySubdomains(osm->n,&osm->iis,&osm->ois)); 9738f3b4b4dSDmitry Karpeev osm->iis = osm->ois = NULL; 9748f3b4b4dSDmitry Karpeev osm->n = n; 9758f3b4b4dSDmitry Karpeev osm->N = PETSC_DETERMINE; 9768f3b4b4dSDmitry Karpeev osm->nmax = PETSC_DETERMINE; 977a35b7b57SDmitry Karpeev if (ois) { 9785f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(n,&osm->ois)); 9798f3b4b4dSDmitry Karpeev for (i=0; i<n; i++) { 9805f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectReference((PetscObject)ois[i])); 9818f3b4b4dSDmitry Karpeev osm->ois[i] = ois[i]; 9828f3b4b4dSDmitry Karpeev } 9838f3b4b4dSDmitry Karpeev /* 9848f3b4b4dSDmitry Karpeev Since the user set the outer subdomains, even if nontrivial overlap was requested via PCGASMSetOverlap(), 9858f3b4b4dSDmitry Karpeev it will be ignored. To avoid confusion later on (e.g., when viewing the PC), the overlap size is set to -1. 9868f3b4b4dSDmitry Karpeev */ 987b1a0a8a3SJed Brown osm->overlap = -1; 988670417b2SFande Kong /* inner subdomains must be provided */ 989*28b400f6SJacob Faibussowitsch PetscCheck(iis,PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"inner indices have to be provided "); 990670417b2SFande Kong }/* end if */ 991a35b7b57SDmitry Karpeev if (iis) { 9925f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(n,&osm->iis)); 9938f3b4b4dSDmitry Karpeev for (i=0; i<n; i++) { 9945f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectReference((PetscObject)iis[i])); 9958f3b4b4dSDmitry Karpeev osm->iis[i] = iis[i]; 9968f3b4b4dSDmitry Karpeev } 997a35b7b57SDmitry Karpeev if (!ois) { 998390e1bf2SBarry Smith osm->ois = NULL; 999670417b2SFande Kong /* if user does not provide outer indices, we will create the corresponding outer indices using osm->overlap =1 in PCSetUp_GASM */ 1000670417b2SFande Kong } 1001670417b2SFande Kong } 100276bd3646SJed Brown if (PetscDefined(USE_DEBUG)) { 1003670417b2SFande Kong PetscInt j,rstart,rend,*covered,lsize; 1004670417b2SFande Kong const PetscInt *indices; 1005670417b2SFande Kong /* check if the inner indices cover and only cover the local portion of the preconditioning matrix */ 10065f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRange(pc->pmat,&rstart,&rend)); 10075f80ce2aSJacob Faibussowitsch CHKERRQ(PetscCalloc1(rend-rstart,&covered)); 1008670417b2SFande Kong /* check if the current processor owns indices from others */ 1009a35b7b57SDmitry Karpeev for (i=0; i<n; i++) { 10105f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetIndices(osm->iis[i],&indices)); 10115f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetLocalSize(osm->iis[i],&lsize)); 1012670417b2SFande Kong for (j=0; j<lsize; j++) { 10132c71b3e2SJacob Faibussowitsch PetscCheckFalse(indices[j]<rstart || indices[j]>=rend,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"inner subdomains can not own an index %d from other processors", indices[j]); 10142c71b3e2SJacob Faibussowitsch else PetscCheckFalse(covered[indices[j]-rstart]==1,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"inner subdomains can not have an overlapping index %d ",indices[j]); 1015670417b2SFande Kong else covered[indices[j]-rstart] = 1; 1016a35b7b57SDmitry Karpeev } 10175f80ce2aSJacob Faibussowitsch CHKERRQ(ISRestoreIndices(osm->iis[i],&indices)); 10188f3b4b4dSDmitry Karpeev } 1019670417b2SFande Kong /* check if we miss any indices */ 1020670417b2SFande Kong for (i=rstart; i<rend; i++) { 10212c71b3e2SJacob Faibussowitsch PetscCheckFalse(!covered[i-rstart],PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"local entity %d was not covered by inner subdomains",i); 1022a35b7b57SDmitry Karpeev } 10235f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(covered)); 1024a35b7b57SDmitry Karpeev } 10258f3b4b4dSDmitry Karpeev if (iis) osm->user_subdomains = PETSC_TRUE; 1026b1a0a8a3SJed Brown PetscFunctionReturn(0); 1027b1a0a8a3SJed Brown } 1028b1a0a8a3SJed Brown 1029f7a08781SBarry Smith static PetscErrorCode PCGASMSetOverlap_GASM(PC pc,PetscInt ovl) 1030b1a0a8a3SJed Brown { 1031f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 1032b1a0a8a3SJed Brown 1033b1a0a8a3SJed Brown PetscFunctionBegin; 10342c71b3e2SJacob Faibussowitsch PetscCheckFalse(ovl < 0,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_OUTOFRANGE,"Negative overlap value requested"); 10352c71b3e2SJacob Faibussowitsch PetscCheckFalse(pc->setupcalled && ovl != osm->overlap,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"PCGASMSetOverlap() should be called before PCSetUp()."); 10362fa5cd67SKarl Rupp if (!pc->setupcalled) osm->overlap = ovl; 1037b1a0a8a3SJed Brown PetscFunctionReturn(0); 1038b1a0a8a3SJed Brown } 1039b1a0a8a3SJed Brown 1040f7a08781SBarry Smith static PetscErrorCode PCGASMSetType_GASM(PC pc,PCGASMType type) 1041b1a0a8a3SJed Brown { 1042f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 1043b1a0a8a3SJed Brown 1044b1a0a8a3SJed Brown PetscFunctionBegin; 1045b1a0a8a3SJed Brown osm->type = type; 1046b1a0a8a3SJed Brown osm->type_set = PETSC_TRUE; 1047b1a0a8a3SJed Brown PetscFunctionReturn(0); 1048b1a0a8a3SJed Brown } 1049b1a0a8a3SJed Brown 1050f7a08781SBarry Smith static PetscErrorCode PCGASMSetSortIndices_GASM(PC pc,PetscBool doSort) 1051b1a0a8a3SJed Brown { 1052f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 1053b1a0a8a3SJed Brown 1054b1a0a8a3SJed Brown PetscFunctionBegin; 1055b1a0a8a3SJed Brown osm->sort_indices = doSort; 1056b1a0a8a3SJed Brown PetscFunctionReturn(0); 1057b1a0a8a3SJed Brown } 1058b1a0a8a3SJed Brown 105944fe18e5SDmitry Karpeev /* 10608f3b4b4dSDmitry Karpeev FIXME: This routine might need to be modified now that multiple ranks per subdomain are allowed. 106144fe18e5SDmitry Karpeev In particular, it would upset the global subdomain number calculation. 106244fe18e5SDmitry Karpeev */ 1063f7a08781SBarry Smith static PetscErrorCode PCGASMGetSubKSP_GASM(PC pc,PetscInt *n,PetscInt *first,KSP **ksp) 1064b1a0a8a3SJed Brown { 1065f746d493SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 1066b1a0a8a3SJed Brown 1067b1a0a8a3SJed Brown PetscFunctionBegin; 10682c71b3e2SJacob Faibussowitsch PetscCheckFalse(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"); 1069b1a0a8a3SJed Brown 10702fa5cd67SKarl Rupp if (n) *n = osm->n; 1071ab3e923fSDmitry Karpeev if (first) { 10725f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Scan(&osm->n,first,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)pc))); 1073ab3e923fSDmitry Karpeev *first -= osm->n; 1074b1a0a8a3SJed Brown } 1075b1a0a8a3SJed Brown if (ksp) { 1076b1a0a8a3SJed Brown /* Assume that local solves are now different; not necessarily 107706b43e7eSDmitry Karpeev true, though! This flag is used only for PCView_GASM() */ 1078b1a0a8a3SJed Brown *ksp = osm->ksp; 10796a4f0f73SDmitry Karpeev osm->same_subdomain_solvers = PETSC_FALSE; 1080b1a0a8a3SJed Brown } 1081b1a0a8a3SJed Brown PetscFunctionReturn(0); 1082ab3e923fSDmitry Karpeev } /* PCGASMGetSubKSP_GASM() */ 1083b1a0a8a3SJed Brown 1084b1a0a8a3SJed Brown /*@C 108506b43e7eSDmitry Karpeev PCGASMSetSubdomains - Sets the subdomains for this processor 108606b43e7eSDmitry Karpeev for the additive Schwarz preconditioner. 1087b1a0a8a3SJed Brown 10888f3b4b4dSDmitry Karpeev Collective on pc 1089b1a0a8a3SJed Brown 1090b1a0a8a3SJed Brown Input Parameters: 10918f3b4b4dSDmitry Karpeev + pc - the preconditioner object 109206b43e7eSDmitry Karpeev . n - the number of subdomains for this processor 10938f3b4b4dSDmitry Karpeev . iis - the index sets that define the inner subdomains (or NULL for PETSc to determine subdomains) 10948f3b4b4dSDmitry 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) 1095b1a0a8a3SJed Brown 1096b1a0a8a3SJed Brown Notes: 109706b43e7eSDmitry Karpeev The IS indices use the parallel, global numbering of the vector entries. 10986a4f0f73SDmitry Karpeev Inner subdomains are those where the correction is applied. 10996a4f0f73SDmitry Karpeev Outer subdomains are those where the residual necessary to obtain the 11006a4f0f73SDmitry Karpeev corrections is obtained (see PCGASMType for the use of inner/outer subdomains). 11016a4f0f73SDmitry Karpeev Both inner and outer subdomains can extend over several processors. 11026a4f0f73SDmitry Karpeev This processor's portion of a subdomain is known as a local subdomain. 11036a4f0f73SDmitry Karpeev 1104670417b2SFande Kong Inner subdomains can not overlap with each other, do not have any entities from remote processors, 1105670417b2SFande Kong and have to cover the entire local subdomain owned by the current processor. The index sets on each 1106670417b2SFande Kong process should be ordered such that the ith local subdomain is connected to the ith remote subdomain 1107670417b2SFande Kong on another MPI process. 1108670417b2SFande Kong 11096a4f0f73SDmitry Karpeev By default the GASM preconditioner uses 1 (local) subdomain per processor. 11106a4f0f73SDmitry Karpeev 1111b1a0a8a3SJed Brown Level: advanced 1112b1a0a8a3SJed Brown 111396dd8680SPierre Jolivet .seealso: PCGASMSetOverlap(), PCGASMGetSubKSP(), 111406b43e7eSDmitry Karpeev PCGASMCreateSubdomains2D(), PCGASMGetSubdomains() 1115b1a0a8a3SJed Brown @*/ 11166a4f0f73SDmitry Karpeev PetscErrorCode PCGASMSetSubdomains(PC pc,PetscInt n,IS iis[],IS ois[]) 1117b1a0a8a3SJed Brown { 11188f3b4b4dSDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 1119b1a0a8a3SJed Brown 1120b1a0a8a3SJed Brown PetscFunctionBegin; 1121b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 11225f80ce2aSJacob Faibussowitsch CHKERRQ(PetscTryMethod(pc,"PCGASMSetSubdomains_C",(PC,PetscInt,IS[],IS[]),(pc,n,iis,ois))); 11238f3b4b4dSDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 1124b1a0a8a3SJed Brown PetscFunctionReturn(0); 1125b1a0a8a3SJed Brown } 1126b1a0a8a3SJed Brown 1127b1a0a8a3SJed Brown /*@ 1128f746d493SDmitry Karpeev PCGASMSetOverlap - Sets the overlap between a pair of subdomains for the 1129b1a0a8a3SJed Brown additive Schwarz preconditioner. Either all or no processors in the 11308f3b4b4dSDmitry Karpeev pc communicator must call this routine. 1131b1a0a8a3SJed Brown 11328f3b4b4dSDmitry Karpeev Logically Collective on pc 1133b1a0a8a3SJed Brown 1134b1a0a8a3SJed Brown Input Parameters: 1135b1a0a8a3SJed Brown + pc - the preconditioner context 11368f3b4b4dSDmitry Karpeev - ovl - the amount of overlap between subdomains (ovl >= 0, default value = 0) 1137b1a0a8a3SJed Brown 1138b1a0a8a3SJed Brown Options Database Key: 113906b43e7eSDmitry Karpeev . -pc_gasm_overlap <overlap> - Sets overlap 1140b1a0a8a3SJed Brown 1141b1a0a8a3SJed Brown Notes: 114206b43e7eSDmitry Karpeev By default the GASM preconditioner uses 1 subdomain per processor. To use 11438f3b4b4dSDmitry Karpeev multiple subdomain per perocessor or "straddling" subdomains that intersect 11448f3b4b4dSDmitry Karpeev multiple processors use PCGASMSetSubdomains() (or option -pc_gasm_total_subdomains <n>). 1145b1a0a8a3SJed Brown 11468f3b4b4dSDmitry Karpeev The overlap defaults to 0, so if one desires that no additional 1147b1a0a8a3SJed Brown overlap be computed beyond what may have been set with a call to 11488f3b4b4dSDmitry Karpeev PCGASMSetSubdomains(), then ovl must be set to be 0. In particular, if one does 11498f3b4b4dSDmitry Karpeev not explicitly set the subdomains in application code, then all overlap would be computed 1150f746d493SDmitry Karpeev internally by PETSc, and using an overlap of 0 would result in an GASM 1151b1a0a8a3SJed Brown variant that is equivalent to the block Jacobi preconditioner. 1152b1a0a8a3SJed Brown 1153b1a0a8a3SJed Brown Note that one can define initial index sets with any overlap via 115406b43e7eSDmitry Karpeev PCGASMSetSubdomains(); the routine PCGASMSetOverlap() merely allows 115506b43e7eSDmitry Karpeev PETSc to extend that overlap further, if desired. 1156b1a0a8a3SJed Brown 1157b1a0a8a3SJed Brown Level: intermediate 1158b1a0a8a3SJed Brown 11598f3b4b4dSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMGetSubKSP(), 116006b43e7eSDmitry Karpeev PCGASMCreateSubdomains2D(), PCGASMGetSubdomains() 1161b1a0a8a3SJed Brown @*/ 11627087cfbeSBarry Smith PetscErrorCode PCGASMSetOverlap(PC pc,PetscInt ovl) 1163b1a0a8a3SJed Brown { 11648f3b4b4dSDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 1165b1a0a8a3SJed Brown 1166b1a0a8a3SJed Brown PetscFunctionBegin; 1167b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1168b1a0a8a3SJed Brown PetscValidLogicalCollectiveInt(pc,ovl,2); 11695f80ce2aSJacob Faibussowitsch CHKERRQ(PetscTryMethod(pc,"PCGASMSetOverlap_C",(PC,PetscInt),(pc,ovl))); 11708f3b4b4dSDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 1171b1a0a8a3SJed Brown PetscFunctionReturn(0); 1172b1a0a8a3SJed Brown } 1173b1a0a8a3SJed Brown 1174b1a0a8a3SJed Brown /*@ 1175f746d493SDmitry Karpeev PCGASMSetType - Sets the type of restriction and interpolation used 1176b1a0a8a3SJed Brown for local problems in the additive Schwarz method. 1177b1a0a8a3SJed Brown 1178b1a0a8a3SJed Brown Logically Collective on PC 1179b1a0a8a3SJed Brown 1180b1a0a8a3SJed Brown Input Parameters: 1181b1a0a8a3SJed Brown + pc - the preconditioner context 1182f746d493SDmitry Karpeev - type - variant of GASM, one of 1183b1a0a8a3SJed Brown .vb 1184f746d493SDmitry Karpeev PC_GASM_BASIC - full interpolation and restriction 1185f746d493SDmitry Karpeev PC_GASM_RESTRICT - full restriction, local processor interpolation 1186f746d493SDmitry Karpeev PC_GASM_INTERPOLATE - full interpolation, local processor restriction 1187f746d493SDmitry Karpeev PC_GASM_NONE - local processor restriction and interpolation 1188b1a0a8a3SJed Brown .ve 1189b1a0a8a3SJed Brown 1190b1a0a8a3SJed Brown Options Database Key: 1191f746d493SDmitry Karpeev . -pc_gasm_type [basic,restrict,interpolate,none] - Sets GASM type 1192b1a0a8a3SJed Brown 1193b1a0a8a3SJed Brown Level: intermediate 1194b1a0a8a3SJed Brown 11958f3b4b4dSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMGetSubKSP(), 1196f746d493SDmitry Karpeev PCGASMCreateSubdomains2D() 1197b1a0a8a3SJed Brown @*/ 11987087cfbeSBarry Smith PetscErrorCode PCGASMSetType(PC pc,PCGASMType type) 1199b1a0a8a3SJed Brown { 1200b1a0a8a3SJed Brown PetscFunctionBegin; 1201b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1202b1a0a8a3SJed Brown PetscValidLogicalCollectiveEnum(pc,type,2); 12035f80ce2aSJacob Faibussowitsch CHKERRQ(PetscTryMethod(pc,"PCGASMSetType_C",(PC,PCGASMType),(pc,type))); 1204b1a0a8a3SJed Brown PetscFunctionReturn(0); 1205b1a0a8a3SJed Brown } 1206b1a0a8a3SJed Brown 1207b1a0a8a3SJed Brown /*@ 1208f746d493SDmitry Karpeev PCGASMSetSortIndices - Determines whether subdomain indices are sorted. 1209b1a0a8a3SJed Brown 1210b1a0a8a3SJed Brown Logically Collective on PC 1211b1a0a8a3SJed Brown 1212b1a0a8a3SJed Brown Input Parameters: 1213b1a0a8a3SJed Brown + pc - the preconditioner context 1214b1a0a8a3SJed Brown - doSort - sort the subdomain indices 1215b1a0a8a3SJed Brown 1216b1a0a8a3SJed Brown Level: intermediate 1217b1a0a8a3SJed Brown 12188f3b4b4dSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMGetSubKSP(), 1219f746d493SDmitry Karpeev PCGASMCreateSubdomains2D() 1220b1a0a8a3SJed Brown @*/ 12217087cfbeSBarry Smith PetscErrorCode PCGASMSetSortIndices(PC pc,PetscBool doSort) 1222b1a0a8a3SJed Brown { 1223b1a0a8a3SJed Brown PetscFunctionBegin; 1224b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1225b1a0a8a3SJed Brown PetscValidLogicalCollectiveBool(pc,doSort,2); 12265f80ce2aSJacob Faibussowitsch CHKERRQ(PetscTryMethod(pc,"PCGASMSetSortIndices_C",(PC,PetscBool),(pc,doSort))); 1227b1a0a8a3SJed Brown PetscFunctionReturn(0); 1228b1a0a8a3SJed Brown } 1229b1a0a8a3SJed Brown 1230b1a0a8a3SJed Brown /*@C 1231f746d493SDmitry Karpeev PCGASMGetSubKSP - Gets the local KSP contexts for all blocks on 1232b1a0a8a3SJed Brown this processor. 1233b1a0a8a3SJed Brown 1234b1a0a8a3SJed Brown Collective on PC iff first_local is requested 1235b1a0a8a3SJed Brown 1236b1a0a8a3SJed Brown Input Parameter: 1237b1a0a8a3SJed Brown . pc - the preconditioner context 1238b1a0a8a3SJed Brown 1239b1a0a8a3SJed Brown Output Parameters: 12400298fd71SBarry Smith + n_local - the number of blocks on this processor or NULL 12410298fd71SBarry Smith . first_local - the global number of the first block on this processor or NULL, 12420298fd71SBarry Smith all processors must request or all must pass NULL 1243b1a0a8a3SJed Brown - ksp - the array of KSP contexts 1244b1a0a8a3SJed Brown 1245b1a0a8a3SJed Brown Note: 1246f746d493SDmitry Karpeev After PCGASMGetSubKSP() the array of KSPes is not to be freed 1247b1a0a8a3SJed Brown 1248b1a0a8a3SJed Brown Currently for some matrix implementations only 1 block per processor 1249b1a0a8a3SJed Brown is supported. 1250b1a0a8a3SJed Brown 1251f746d493SDmitry Karpeev You must call KSPSetUp() before calling PCGASMGetSubKSP(). 1252b1a0a8a3SJed Brown 1253b1a0a8a3SJed Brown Level: advanced 1254b1a0a8a3SJed Brown 12558f3b4b4dSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMSetOverlap(), 1256f746d493SDmitry Karpeev PCGASMCreateSubdomains2D(), 1257b1a0a8a3SJed Brown @*/ 12587087cfbeSBarry Smith PetscErrorCode PCGASMGetSubKSP(PC pc,PetscInt *n_local,PetscInt *first_local,KSP *ksp[]) 1259b1a0a8a3SJed Brown { 1260b1a0a8a3SJed Brown PetscFunctionBegin; 1261b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 12625f80ce2aSJacob Faibussowitsch CHKERRQ(PetscUseMethod(pc,"PCGASMGetSubKSP_C",(PC,PetscInt*,PetscInt*,KSP **),(pc,n_local,first_local,ksp))); 1263b1a0a8a3SJed Brown PetscFunctionReturn(0); 1264b1a0a8a3SJed Brown } 1265b1a0a8a3SJed Brown 1266b1a0a8a3SJed Brown /* -------------------------------------------------------------------------------------*/ 1267b1a0a8a3SJed Brown /*MC 1268f746d493SDmitry Karpeev PCGASM - Use the (restricted) additive Schwarz method, each block is (approximately) solved with 1269b1a0a8a3SJed Brown its own KSP object. 1270b1a0a8a3SJed Brown 1271b1a0a8a3SJed Brown Options Database Keys: 12728f3b4b4dSDmitry Karpeev + -pc_gasm_total_subdomains <n> - Sets total number of local subdomains to be distributed among processors 127306b43e7eSDmitry Karpeev . -pc_gasm_view_subdomains - activates the printing of subdomain indices in PCView(), -ksp_view or -snes_view 127406b43e7eSDmitry Karpeev . -pc_gasm_print_subdomains - activates the printing of subdomain indices in PCSetUp() 127506b43e7eSDmitry Karpeev . -pc_gasm_overlap <ovl> - Sets overlap by which to (automatically) extend local subdomains 1276f746d493SDmitry Karpeev - -pc_gasm_type [basic,restrict,interpolate,none] - Sets GASM type 1277b1a0a8a3SJed Brown 1278b1a0a8a3SJed Brown IMPORTANT: If you run with, for example, 3 blocks on 1 processor or 3 blocks on 3 processors you 1279f746d493SDmitry Karpeev will get a different convergence rate due to the default option of -pc_gasm_type restrict. Use 1280f746d493SDmitry Karpeev -pc_gasm_type basic to use the standard GASM. 1281b1a0a8a3SJed Brown 128295452b02SPatrick Sanan Notes: 128395452b02SPatrick Sanan Blocks can be shared by multiple processes. 1284b1a0a8a3SJed Brown 1285b1a0a8a3SJed Brown To set options on the solvers for each block append -sub_ to all the KSP, and PC 1286b1a0a8a3SJed Brown options database keys. For example, -sub_pc_type ilu -sub_pc_factor_levels 1 -sub_ksp_type preonly 1287b1a0a8a3SJed Brown 1288f746d493SDmitry Karpeev To set the options on the solvers separate for each block call PCGASMGetSubKSP() 1289b1a0a8a3SJed Brown and set the options directly on the resulting KSP object (you can access its PC 1290b1a0a8a3SJed Brown with KSPGetPC()) 1291b1a0a8a3SJed Brown 1292b1a0a8a3SJed Brown Level: beginner 1293b1a0a8a3SJed Brown 1294b1a0a8a3SJed Brown References: 1295606c0280SSatish Balay + * - M Dryja, OB Widlund, An additive variant of the Schwarz alternating method for the case of many subregions 129696a0c994SBarry Smith Courant Institute, New York University Technical report 1297606c0280SSatish Balay - * - Barry Smith, Petter Bjorstad, and William Gropp, Domain Decompositions: Parallel Multilevel Methods for Elliptic Partial Differential Equations, 129896a0c994SBarry Smith Cambridge University Press. 1299b1a0a8a3SJed Brown 1300b1a0a8a3SJed Brown .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, 130149517cdeSBarry Smith PCBJACOBI, PCGASMGetSubKSP(), PCGASMSetSubdomains(), 130234a84908Sprj- PCSetModifySubMatrices(), PCGASMSetOverlap(), PCGASMSetType() 1303b1a0a8a3SJed Brown 1304b1a0a8a3SJed Brown M*/ 1305b1a0a8a3SJed Brown 13068cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_GASM(PC pc) 1307b1a0a8a3SJed Brown { 1308f746d493SDmitry Karpeev PC_GASM *osm; 1309b1a0a8a3SJed Brown 1310b1a0a8a3SJed Brown PetscFunctionBegin; 13115f80ce2aSJacob Faibussowitsch CHKERRQ(PetscNewLog(pc,&osm)); 13122fa5cd67SKarl Rupp 13138f3b4b4dSDmitry Karpeev osm->N = PETSC_DETERMINE; 131406b43e7eSDmitry Karpeev osm->n = PETSC_DECIDE; 13158f3b4b4dSDmitry Karpeev osm->nmax = PETSC_DETERMINE; 13168f3b4b4dSDmitry Karpeev osm->overlap = 0; 13170a545947SLisandro Dalcin osm->ksp = NULL; 13180a545947SLisandro Dalcin osm->gorestriction = NULL; 13190a545947SLisandro Dalcin osm->girestriction = NULL; 13200a545947SLisandro Dalcin osm->pctoouter = NULL; 13210a545947SLisandro Dalcin osm->gx = NULL; 13220a545947SLisandro Dalcin osm->gy = NULL; 13230a545947SLisandro Dalcin osm->x = NULL; 13240a545947SLisandro Dalcin osm->y = NULL; 13250a545947SLisandro Dalcin osm->pcx = NULL; 13260a545947SLisandro Dalcin osm->pcy = NULL; 13270a545947SLisandro Dalcin osm->permutationIS = NULL; 13280a545947SLisandro Dalcin osm->permutationP = NULL; 13290a545947SLisandro Dalcin osm->pcmat = NULL; 13300a545947SLisandro Dalcin osm->ois = NULL; 13310a545947SLisandro Dalcin osm->iis = NULL; 13320a545947SLisandro Dalcin osm->pmat = NULL; 1333f746d493SDmitry Karpeev osm->type = PC_GASM_RESTRICT; 13346a4f0f73SDmitry Karpeev osm->same_subdomain_solvers = PETSC_TRUE; 1335b1a0a8a3SJed Brown osm->sort_indices = PETSC_TRUE; 1336d709ea83SDmitry Karpeev osm->dm_subdomains = PETSC_FALSE; 1337ea91fabdSFande Kong osm->hierarchicalpartitioning = PETSC_FALSE; 1338b1a0a8a3SJed Brown 1339b1a0a8a3SJed Brown pc->data = (void*)osm; 1340f746d493SDmitry Karpeev pc->ops->apply = PCApply_GASM; 134148e38a8aSPierre Jolivet pc->ops->matapply = PCMatApply_GASM; 1342f746d493SDmitry Karpeev pc->ops->applytranspose = PCApplyTranspose_GASM; 1343f746d493SDmitry Karpeev pc->ops->setup = PCSetUp_GASM; 1344a06653b4SBarry Smith pc->ops->reset = PCReset_GASM; 1345f746d493SDmitry Karpeev pc->ops->destroy = PCDestroy_GASM; 1346f746d493SDmitry Karpeev pc->ops->setfromoptions = PCSetFromOptions_GASM; 1347f746d493SDmitry Karpeev pc->ops->setuponblocks = PCSetUpOnBlocks_GASM; 1348f746d493SDmitry Karpeev pc->ops->view = PCView_GASM; 13490a545947SLisandro Dalcin pc->ops->applyrichardson = NULL; 1350b1a0a8a3SJed Brown 13515f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetSubdomains_C",PCGASMSetSubdomains_GASM)); 13525f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetOverlap_C",PCGASMSetOverlap_GASM)); 13535f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetType_C",PCGASMSetType_GASM)); 13545f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)pc,"PCGASMSetSortIndices_C",PCGASMSetSortIndices_GASM)); 13555f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectComposeFunction((PetscObject)pc,"PCGASMGetSubKSP_C",PCGASMGetSubKSP_GASM)); 1356b1a0a8a3SJed Brown PetscFunctionReturn(0); 1357b1a0a8a3SJed Brown } 1358b1a0a8a3SJed Brown 13598f3b4b4dSDmitry Karpeev PetscErrorCode PCGASMCreateLocalSubdomains(Mat A, PetscInt nloc, IS *iis[]) 1360b1a0a8a3SJed Brown { 1361b1a0a8a3SJed Brown MatPartitioning mpart; 1362b1a0a8a3SJed Brown const char *prefix; 1363b1a0a8a3SJed Brown PetscInt i,j,rstart,rend,bs; 1364976e8c5aSLisandro Dalcin PetscBool hasop, isbaij = PETSC_FALSE,foundpart = PETSC_FALSE; 13650298fd71SBarry Smith Mat Ad = NULL, adj; 1366b1a0a8a3SJed Brown IS ispart,isnumb,*is; 1367b1a0a8a3SJed Brown 1368b1a0a8a3SJed Brown PetscFunctionBegin; 13692c71b3e2SJacob Faibussowitsch PetscCheckFalse(nloc < 1,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"number of local subdomains must > 0, got nloc = %D",nloc); 1370b1a0a8a3SJed Brown 1371b1a0a8a3SJed Brown /* Get prefix, row distribution, and block size */ 13725f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOptionsPrefix(A,&prefix)); 13735f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRange(A,&rstart,&rend)); 13745f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetBlockSize(A,&bs)); 13752c71b3e2SJacob Faibussowitsch PetscCheckFalse(rstart/bs*bs != rstart || rend/bs*bs != rend,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"bad row distribution [%D,%D) for matrix block size %D",rstart,rend,bs); 1376b1a0a8a3SJed Brown 1377b1a0a8a3SJed Brown /* Get diagonal block from matrix if possible */ 13785f80ce2aSJacob Faibussowitsch CHKERRQ(MatHasOperation(A,MATOP_GET_DIAGONAL_BLOCK,&hasop)); 1379976e8c5aSLisandro Dalcin if (hasop) { 13805f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetDiagonalBlock(A,&Ad)); 1381b1a0a8a3SJed Brown } 1382b1a0a8a3SJed Brown if (Ad) { 13835f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectBaseTypeCompare((PetscObject)Ad,MATSEQBAIJ,&isbaij)); 13845f80ce2aSJacob Faibussowitsch if (!isbaij) CHKERRQ(PetscObjectBaseTypeCompare((PetscObject)Ad,MATSEQSBAIJ,&isbaij)); 1385b1a0a8a3SJed Brown } 13868f3b4b4dSDmitry Karpeev if (Ad && nloc > 1) { 1387b1a0a8a3SJed Brown PetscBool match,done; 1388b1a0a8a3SJed Brown /* Try to setup a good matrix partitioning if available */ 13895f80ce2aSJacob Faibussowitsch CHKERRQ(MatPartitioningCreate(PETSC_COMM_SELF,&mpart)); 13905f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectSetOptionsPrefix((PetscObject)mpart,prefix)); 13915f80ce2aSJacob Faibussowitsch CHKERRQ(MatPartitioningSetFromOptions(mpart)); 13925f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGCURRENT,&match)); 1393b1a0a8a3SJed Brown if (!match) { 13945f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectTypeCompare((PetscObject)mpart,MATPARTITIONINGSQUARE,&match)); 1395b1a0a8a3SJed Brown } 1396b1a0a8a3SJed Brown if (!match) { /* assume a "good" partitioner is available */ 13971a83f524SJed Brown PetscInt na; 13981a83f524SJed Brown const PetscInt *ia,*ja; 13995f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done)); 1400b1a0a8a3SJed Brown if (done) { 1401b1a0a8a3SJed Brown /* Build adjacency matrix by hand. Unfortunately a call to 1402b1a0a8a3SJed Brown MatConvert(Ad,MATMPIADJ,MAT_INITIAL_MATRIX,&adj) will 1403b1a0a8a3SJed Brown remove the block-aij structure and we cannot expect 1404b1a0a8a3SJed Brown MatPartitioning to split vertices as we need */ 14050a545947SLisandro Dalcin PetscInt i,j,len,nnz,cnt,*iia=NULL,*jja=NULL; 14061a83f524SJed Brown const PetscInt *row; 1407b1a0a8a3SJed Brown nnz = 0; 1408b1a0a8a3SJed Brown for (i=0; i<na; i++) { /* count number of nonzeros */ 1409b1a0a8a3SJed Brown len = ia[i+1] - ia[i]; 1410b1a0a8a3SJed Brown row = ja + ia[i]; 1411b1a0a8a3SJed Brown for (j=0; j<len; j++) { 1412b1a0a8a3SJed Brown if (row[j] == i) { /* don't count diagonal */ 1413b1a0a8a3SJed Brown len--; break; 1414b1a0a8a3SJed Brown } 1415b1a0a8a3SJed Brown } 1416b1a0a8a3SJed Brown nnz += len; 1417b1a0a8a3SJed Brown } 14185f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(na+1,&iia)); 14195f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(nnz,&jja)); 1420b1a0a8a3SJed Brown nnz = 0; 1421b1a0a8a3SJed Brown iia[0] = 0; 1422b1a0a8a3SJed Brown for (i=0; i<na; i++) { /* fill adjacency */ 1423b1a0a8a3SJed Brown cnt = 0; 1424b1a0a8a3SJed Brown len = ia[i+1] - ia[i]; 1425b1a0a8a3SJed Brown row = ja + ia[i]; 1426b1a0a8a3SJed Brown for (j=0; j<len; j++) { 14272fa5cd67SKarl Rupp if (row[j] != i) jja[nnz+cnt++] = row[j]; /* if not diagonal */ 1428b1a0a8a3SJed Brown } 1429b1a0a8a3SJed Brown nnz += cnt; 1430b1a0a8a3SJed Brown iia[i+1] = nnz; 1431b1a0a8a3SJed Brown } 1432b1a0a8a3SJed Brown /* Partitioning of the adjacency matrix */ 14335f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateMPIAdj(PETSC_COMM_SELF,na,na,iia,jja,NULL,&adj)); 14345f80ce2aSJacob Faibussowitsch CHKERRQ(MatPartitioningSetAdjacency(mpart,adj)); 14355f80ce2aSJacob Faibussowitsch CHKERRQ(MatPartitioningSetNParts(mpart,nloc)); 14365f80ce2aSJacob Faibussowitsch CHKERRQ(MatPartitioningApply(mpart,&ispart)); 14375f80ce2aSJacob Faibussowitsch CHKERRQ(ISPartitioningToNumbering(ispart,&isnumb)); 14385f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&adj)); 1439b1a0a8a3SJed Brown foundpart = PETSC_TRUE; 1440b1a0a8a3SJed Brown } 14415f80ce2aSJacob Faibussowitsch CHKERRQ(MatRestoreRowIJ(Ad,0,PETSC_TRUE,isbaij,&na,&ia,&ja,&done)); 1442b1a0a8a3SJed Brown } 14435f80ce2aSJacob Faibussowitsch CHKERRQ(MatPartitioningDestroy(&mpart)); 1444b1a0a8a3SJed Brown } 14455f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(nloc,&is)); 1446b1a0a8a3SJed Brown if (!foundpart) { 1447b1a0a8a3SJed Brown 1448b1a0a8a3SJed Brown /* Partitioning by contiguous chunks of rows */ 1449b1a0a8a3SJed Brown 1450b1a0a8a3SJed Brown PetscInt mbs = (rend-rstart)/bs; 1451b1a0a8a3SJed Brown PetscInt start = rstart; 14528f3b4b4dSDmitry Karpeev for (i=0; i<nloc; i++) { 14538f3b4b4dSDmitry Karpeev PetscInt count = (mbs/nloc + ((mbs % nloc) > i)) * bs; 14545f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateStride(PETSC_COMM_SELF,count,start,1,&is[i])); 1455b1a0a8a3SJed Brown start += count; 1456b1a0a8a3SJed Brown } 1457b1a0a8a3SJed Brown 1458b1a0a8a3SJed Brown } else { 1459b1a0a8a3SJed Brown 1460b1a0a8a3SJed Brown /* Partitioning by adjacency of diagonal block */ 1461b1a0a8a3SJed Brown 1462b1a0a8a3SJed Brown const PetscInt *numbering; 1463b1a0a8a3SJed Brown PetscInt *count,nidx,*indices,*newidx,start=0; 1464b1a0a8a3SJed Brown /* Get node count in each partition */ 14655f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(nloc,&count)); 14665f80ce2aSJacob Faibussowitsch CHKERRQ(ISPartitioningCount(ispart,nloc,count)); 1467b1a0a8a3SJed Brown if (isbaij && bs > 1) { /* adjust for the block-aij case */ 14688f3b4b4dSDmitry Karpeev for (i=0; i<nloc; i++) count[i] *= bs; 1469b1a0a8a3SJed Brown } 1470b1a0a8a3SJed Brown /* Build indices from node numbering */ 14715f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetLocalSize(isnumb,&nidx)); 14725f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(nidx,&indices)); 1473b1a0a8a3SJed Brown for (i=0; i<nidx; i++) indices[i] = i; /* needs to be initialized */ 14745f80ce2aSJacob Faibussowitsch CHKERRQ(ISGetIndices(isnumb,&numbering)); 14755f80ce2aSJacob Faibussowitsch CHKERRQ(PetscSortIntWithPermutation(nidx,numbering,indices)); 14765f80ce2aSJacob Faibussowitsch CHKERRQ(ISRestoreIndices(isnumb,&numbering)); 1477b1a0a8a3SJed Brown if (isbaij && bs > 1) { /* adjust for the block-aij case */ 14785f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(nidx*bs,&newidx)); 14792fa5cd67SKarl Rupp for (i=0; i<nidx; i++) { 14802fa5cd67SKarl Rupp for (j=0; j<bs; j++) newidx[i*bs+j] = indices[i]*bs + j; 14812fa5cd67SKarl Rupp } 14825f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(indices)); 1483b1a0a8a3SJed Brown nidx *= bs; 1484b1a0a8a3SJed Brown indices = newidx; 1485b1a0a8a3SJed Brown } 1486b1a0a8a3SJed Brown /* Shift to get global indices */ 1487b1a0a8a3SJed Brown for (i=0; i<nidx; i++) indices[i] += rstart; 1488b1a0a8a3SJed Brown 1489b1a0a8a3SJed Brown /* Build the index sets for each block */ 14908f3b4b4dSDmitry Karpeev for (i=0; i<nloc; i++) { 14915f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateGeneral(PETSC_COMM_SELF,count[i],&indices[start],PETSC_COPY_VALUES,&is[i])); 14925f80ce2aSJacob Faibussowitsch CHKERRQ(ISSort(is[i])); 1493b1a0a8a3SJed Brown start += count[i]; 1494b1a0a8a3SJed Brown } 1495b1a0a8a3SJed Brown 14965f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(count)); 14975f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree(indices)); 14985f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&isnumb)); 14995f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&ispart)); 1500b1a0a8a3SJed Brown } 15016a4f0f73SDmitry Karpeev *iis = is; 15028f3b4b4dSDmitry Karpeev PetscFunctionReturn(0); 15038f3b4b4dSDmitry Karpeev } 15048f3b4b4dSDmitry Karpeev 1505b541e6a4SDmitry Karpeev PETSC_INTERN PetscErrorCode PCGASMCreateStraddlingSubdomains(Mat A,PetscInt N,PetscInt *n,IS *iis[]) 15068f3b4b4dSDmitry Karpeev { 15078f3b4b4dSDmitry Karpeev PetscFunctionBegin; 15085f80ce2aSJacob Faibussowitsch CHKERRQ(MatSubdomainsCreateCoalesce(A,N,n,iis)); 15098f3b4b4dSDmitry Karpeev PetscFunctionReturn(0); 15108f3b4b4dSDmitry Karpeev } 15118f3b4b4dSDmitry Karpeev 15128f3b4b4dSDmitry Karpeev /*@C 15138f3b4b4dSDmitry Karpeev PCGASMCreateSubdomains - Creates n index sets defining n nonoverlapping subdomains for the additive 15148f3b4b4dSDmitry Karpeev Schwarz preconditioner for a any problem based on its matrix. 15158f3b4b4dSDmitry Karpeev 15168f3b4b4dSDmitry Karpeev Collective 15178f3b4b4dSDmitry Karpeev 15188f3b4b4dSDmitry Karpeev Input Parameters: 15198f3b4b4dSDmitry Karpeev + A - The global matrix operator 15208f3b4b4dSDmitry Karpeev - N - the number of global subdomains requested 15218f3b4b4dSDmitry Karpeev 15228f3b4b4dSDmitry Karpeev Output Parameters: 15238f3b4b4dSDmitry Karpeev + n - the number of subdomains created on this processor 15248f3b4b4dSDmitry Karpeev - iis - the array of index sets defining the local inner subdomains (on which the correction is applied) 15258f3b4b4dSDmitry Karpeev 15268f3b4b4dSDmitry Karpeev Level: advanced 15278f3b4b4dSDmitry Karpeev 15288f3b4b4dSDmitry Karpeev Note: When N >= A's communicator size, each subdomain is local -- contained within a single processor. 15298f3b4b4dSDmitry Karpeev When N < size, the subdomains are 'straddling' (processor boundaries) and are no longer local. 15308f3b4b4dSDmitry Karpeev The resulting subdomains can be use in PCGASMSetSubdomains(pc,n,iss,NULL). The overlapping 15318f3b4b4dSDmitry Karpeev outer subdomains will be automatically generated from these according to the requested amount of 15328f3b4b4dSDmitry Karpeev overlap; this is currently supported only with local subdomains. 15338f3b4b4dSDmitry Karpeev 15348f3b4b4dSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMDestroySubdomains() 15358f3b4b4dSDmitry Karpeev @*/ 1536b541e6a4SDmitry Karpeev PetscErrorCode PCGASMCreateSubdomains(Mat A,PetscInt N,PetscInt *n,IS *iis[]) 15378f3b4b4dSDmitry Karpeev { 15388f3b4b4dSDmitry Karpeev PetscMPIInt size; 15398f3b4b4dSDmitry Karpeev 15408f3b4b4dSDmitry Karpeev PetscFunctionBegin; 15418f3b4b4dSDmitry Karpeev PetscValidHeaderSpecific(A,MAT_CLASSID,1); 15428f3b4b4dSDmitry Karpeev PetscValidPointer(iis,4); 15438f3b4b4dSDmitry Karpeev 15442c71b3e2SJacob Faibussowitsch PetscCheckFalse(N < 1,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Number of subdomains must be > 0, N = %D",N); 15455f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PetscObjectComm((PetscObject)A),&size)); 15468f3b4b4dSDmitry Karpeev if (N >= size) { 15478f3b4b4dSDmitry Karpeev *n = N/size + (N%size); 15485f80ce2aSJacob Faibussowitsch CHKERRQ(PCGASMCreateLocalSubdomains(A,*n,iis)); 15496a4f0f73SDmitry Karpeev } else { 15505f80ce2aSJacob Faibussowitsch CHKERRQ(PCGASMCreateStraddlingSubdomains(A,N,n,iis)); 15516a4f0f73SDmitry Karpeev } 1552b1a0a8a3SJed Brown PetscFunctionReturn(0); 1553b1a0a8a3SJed Brown } 1554b1a0a8a3SJed Brown 1555b1a0a8a3SJed Brown /*@C 1556f746d493SDmitry Karpeev PCGASMDestroySubdomains - Destroys the index sets created with 15578f3b4b4dSDmitry Karpeev PCGASMCreateSubdomains() or PCGASMCreateSubdomains2D. Should be 155806b43e7eSDmitry Karpeev called after setting subdomains with PCGASMSetSubdomains(). 1559b1a0a8a3SJed Brown 1560b1a0a8a3SJed Brown Collective 1561b1a0a8a3SJed Brown 1562b1a0a8a3SJed Brown Input Parameters: 1563b1a0a8a3SJed Brown + n - the number of index sets 15646a4f0f73SDmitry Karpeev . iis - the array of inner subdomains, 15650298fd71SBarry Smith - ois - the array of outer subdomains, can be NULL 1566b1a0a8a3SJed Brown 15676a4f0f73SDmitry Karpeev Level: intermediate 15686a4f0f73SDmitry Karpeev 156995452b02SPatrick Sanan Notes: 157095452b02SPatrick Sanan this is merely a convenience subroutine that walks each list, 15712c112581SDmitry Karpeev destroys each IS on the list, and then frees the list. At the end the 15722c112581SDmitry Karpeev list pointers are set to NULL. 1573b1a0a8a3SJed Brown 15748f3b4b4dSDmitry Karpeev .seealso: PCGASMCreateSubdomains(), PCGASMSetSubdomains() 1575b1a0a8a3SJed Brown @*/ 15762c112581SDmitry Karpeev PetscErrorCode PCGASMDestroySubdomains(PetscInt n,IS **iis,IS **ois) 1577b1a0a8a3SJed Brown { 1578b1a0a8a3SJed Brown PetscInt i; 15795fd66863SKarl Rupp 1580b1a0a8a3SJed Brown PetscFunctionBegin; 1581a35b7b57SDmitry Karpeev if (n <= 0) PetscFunctionReturn(0); 15826a4f0f73SDmitry Karpeev if (ois) { 15832c112581SDmitry Karpeev PetscValidPointer(ois,3); 15842c112581SDmitry Karpeev if (*ois) { 15852c112581SDmitry Karpeev PetscValidPointer(*ois,3); 1586a35b7b57SDmitry Karpeev for (i=0; i<n; i++) { 15875f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&(*ois)[i])); 1588a35b7b57SDmitry Karpeev } 15895f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree((*ois))); 15902c112581SDmitry Karpeev } 1591b1a0a8a3SJed Brown } 1592b9d0fdaaSFande Kong if (iis) { 1593b9d0fdaaSFande Kong PetscValidPointer(iis,2); 1594b9d0fdaaSFande Kong if (*iis) { 1595b9d0fdaaSFande Kong PetscValidPointer(*iis,2); 1596b9d0fdaaSFande Kong for (i=0; i<n; i++) { 15975f80ce2aSJacob Faibussowitsch CHKERRQ(ISDestroy(&(*iis)[i])); 1598b9d0fdaaSFande Kong } 15995f80ce2aSJacob Faibussowitsch CHKERRQ(PetscFree((*iis))); 1600b9d0fdaaSFande Kong } 1601b9d0fdaaSFande Kong } 1602b1a0a8a3SJed Brown PetscFunctionReturn(0); 1603b1a0a8a3SJed Brown } 1604b1a0a8a3SJed Brown 1605af538404SDmitry Karpeev #define PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,xleft_loc,ylow_loc,xright_loc,yhigh_loc,n) \ 1606af538404SDmitry Karpeev { \ 1607af538404SDmitry Karpeev PetscInt first_row = first/M, last_row = last/M+1; \ 1608af538404SDmitry Karpeev /* \ 1609af538404SDmitry Karpeev Compute ylow_loc and yhigh_loc so that (ylow_loc,xleft) and (yhigh_loc,xright) are the corners \ 1610af538404SDmitry Karpeev of the bounding box of the intersection of the subdomain with the local ownership range (local \ 1611af538404SDmitry Karpeev subdomain). \ 1612af538404SDmitry Karpeev Also compute xleft_loc and xright_loc as the lower and upper bounds on the first and last rows \ 1613af538404SDmitry Karpeev of the intersection. \ 1614af538404SDmitry Karpeev */ \ 1615af538404SDmitry Karpeev /* ylow_loc is the grid row containing the first element of the local sumbdomain */ \ 1616eec7e3faSDmitry Karpeev *ylow_loc = PetscMax(first_row,ylow); \ 1617af538404SDmitry Karpeev /* xleft_loc is the offset of first element of the local subdomain within its grid row (might actually be outside the local subdomain) */ \ 1618eec7e3faSDmitry Karpeev *xleft_loc = *ylow_loc==first_row ? PetscMax(first%M,xleft) : xleft; \ 1619af538404SDmitry Karpeev /* yhigh_loc is the grid row above the last local subdomain element */ \ 1620eec7e3faSDmitry Karpeev *yhigh_loc = PetscMin(last_row,yhigh); \ 1621af538404SDmitry Karpeev /* xright is the offset of the end of the local subdomain within its grid row (might actually be outside the local subdomain) */ \ 1622eec7e3faSDmitry Karpeev *xright_loc = *yhigh_loc==last_row ? PetscMin(xright,last%M) : xright; \ 1623af538404SDmitry Karpeev /* Now compute the size of the local subdomain n. */ \ 1624c3518bceSDmitry Karpeev *n = 0; \ 1625eec7e3faSDmitry Karpeev if (*ylow_loc < *yhigh_loc) { \ 1626af538404SDmitry Karpeev PetscInt width = xright-xleft; \ 1627eec7e3faSDmitry Karpeev *n += width*(*yhigh_loc-*ylow_loc-1); \ 1628eec7e3faSDmitry Karpeev *n += PetscMin(PetscMax(*xright_loc-xleft,0),width); \ 1629eec7e3faSDmitry Karpeev *n -= PetscMin(PetscMax(*xleft_loc-xleft,0), width); \ 1630af538404SDmitry Karpeev } \ 1631af538404SDmitry Karpeev } 1632af538404SDmitry Karpeev 1633b1a0a8a3SJed Brown /*@ 1634f746d493SDmitry Karpeev PCGASMCreateSubdomains2D - Creates the index sets for the overlapping Schwarz 1635b1a0a8a3SJed Brown preconditioner for a two-dimensional problem on a regular grid. 1636b1a0a8a3SJed Brown 1637af538404SDmitry Karpeev Collective 1638b1a0a8a3SJed Brown 1639b1a0a8a3SJed Brown Input Parameters: 16406b867d5aSJose E. Roman + pc - the preconditioner context 16416b867d5aSJose E. Roman . M - the global number of grid points in the x direction 16426b867d5aSJose E. Roman . N - the global number of grid points in the y direction 16436b867d5aSJose E. Roman . Mdomains - the global number of subdomains in the x direction 16446b867d5aSJose E. Roman . Ndomains - the global number of subdomains in the y direction 1645b1a0a8a3SJed Brown . dof - degrees of freedom per node 1646b1a0a8a3SJed Brown - overlap - overlap in mesh lines 1647b1a0a8a3SJed Brown 1648b1a0a8a3SJed Brown Output Parameters: 1649af538404SDmitry Karpeev + Nsub - the number of local subdomains created 16506a4f0f73SDmitry Karpeev . iis - array of index sets defining inner (nonoverlapping) subdomains 16516a4f0f73SDmitry Karpeev - ois - array of index sets defining outer (overlapping, if overlap > 0) subdomains 1652b1a0a8a3SJed Brown 1653b1a0a8a3SJed Brown Level: advanced 1654b1a0a8a3SJed Brown 16558f3b4b4dSDmitry Karpeev .seealso: PCGASMSetSubdomains(), PCGASMGetSubKSP(), PCGASMSetOverlap() 1656b1a0a8a3SJed Brown @*/ 16576a4f0f73SDmitry Karpeev PetscErrorCode PCGASMCreateSubdomains2D(PC pc,PetscInt M,PetscInt N,PetscInt Mdomains,PetscInt Ndomains,PetscInt dof,PetscInt overlap,PetscInt *nsub,IS **iis,IS **ois) 1658b1a0a8a3SJed Brown { 16592bbb417fSDmitry Karpeev PetscMPIInt size, rank; 16602bbb417fSDmitry Karpeev PetscInt i, j; 16612bbb417fSDmitry Karpeev PetscInt maxheight, maxwidth; 16622bbb417fSDmitry Karpeev PetscInt xstart, xleft, xright, xleft_loc, xright_loc; 16632bbb417fSDmitry Karpeev PetscInt ystart, ylow, yhigh, ylow_loc, yhigh_loc; 1664eec7e3faSDmitry Karpeev PetscInt x[2][2], y[2][2], n[2]; 16652bbb417fSDmitry Karpeev PetscInt first, last; 16662bbb417fSDmitry Karpeev PetscInt nidx, *idx; 16672bbb417fSDmitry Karpeev PetscInt ii,jj,s,q,d; 1668af538404SDmitry Karpeev PetscInt k,kk; 16692bbb417fSDmitry Karpeev PetscMPIInt color; 16702bbb417fSDmitry Karpeev MPI_Comm comm, subcomm; 16710a545947SLisandro Dalcin IS **xis = NULL, **is = ois, **is_local = iis; 1672d34fcf5fSBarry Smith 1673b1a0a8a3SJed Brown PetscFunctionBegin; 16745f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectGetComm((PetscObject)pc, &comm)); 16755f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(comm, &size)); 16765f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_rank(comm, &rank)); 16775f80ce2aSJacob Faibussowitsch CHKERRQ(MatGetOwnershipRange(pc->pmat, &first, &last)); 16782c71b3e2SJacob Faibussowitsch PetscCheckFalse(first%dof || last%dof,PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Matrix row partitioning unsuitable for domain decomposition: local row range (%D,%D) " 16792bbb417fSDmitry Karpeev "does not respect the number of degrees of freedom per grid point %D", first, last, dof); 1680d34fcf5fSBarry Smith 1681af538404SDmitry Karpeev /* Determine the number of domains with nonzero intersections with the local ownership range. */ 16822bbb417fSDmitry Karpeev s = 0; 1683af538404SDmitry Karpeev ystart = 0; 1684af538404SDmitry Karpeev for (j=0; j<Ndomains; ++j) { 1685af538404SDmitry Karpeev maxheight = N/Ndomains + ((N % Ndomains) > j); /* Maximal height of subdomain */ 16862c71b3e2SJacob Faibussowitsch PetscCheckFalse(maxheight < 2,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many %D subdomains in the vertical directon for mesh height %D", Ndomains, N); 1687eec7e3faSDmitry Karpeev /* Vertical domain limits with an overlap. */ 1688eec7e3faSDmitry Karpeev ylow = PetscMax(ystart - overlap,0); 1689eec7e3faSDmitry Karpeev yhigh = PetscMin(ystart + maxheight + overlap,N); 1690b1a0a8a3SJed Brown xstart = 0; 1691af538404SDmitry Karpeev for (i=0; i<Mdomains; ++i) { 1692af538404SDmitry Karpeev maxwidth = M/Mdomains + ((M % Mdomains) > i); /* Maximal width of subdomain */ 16932c71b3e2SJacob Faibussowitsch PetscCheckFalse(maxwidth < 2,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many %D subdomains in the horizontal direction for mesh width %D", Mdomains, M); 1694eec7e3faSDmitry Karpeev /* Horizontal domain limits with an overlap. */ 1695eec7e3faSDmitry Karpeev xleft = PetscMax(xstart - overlap,0); 1696eec7e3faSDmitry Karpeev xright = PetscMin(xstart + maxwidth + overlap,M); 1697af538404SDmitry Karpeev /* 1698af538404SDmitry Karpeev Determine whether this subdomain intersects this processor's ownership range of pc->pmat. 1699af538404SDmitry Karpeev */ 1700c3518bceSDmitry Karpeev PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,(&xleft_loc),(&ylow_loc),(&xright_loc),(&yhigh_loc),(&nidx)); 17012fa5cd67SKarl Rupp if (nidx) ++s; 1702af538404SDmitry Karpeev xstart += maxwidth; 1703af538404SDmitry Karpeev } /* for (i = 0; i < Mdomains; ++i) */ 1704af538404SDmitry Karpeev ystart += maxheight; 1705af538404SDmitry Karpeev } /* for (j = 0; j < Ndomains; ++j) */ 17062fa5cd67SKarl Rupp 1707af538404SDmitry Karpeev /* Now we can allocate the necessary number of ISs. */ 1708af538404SDmitry Karpeev *nsub = s; 17095f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(*nsub,is)); 17105f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(*nsub,is_local)); 1711af538404SDmitry Karpeev s = 0; 1712af538404SDmitry Karpeev ystart = 0; 1713af538404SDmitry Karpeev for (j=0; j<Ndomains; ++j) { 1714af538404SDmitry Karpeev maxheight = N/Ndomains + ((N % Ndomains) > j); /* Maximal height of subdomain */ 17152c71b3e2SJacob Faibussowitsch PetscCheckFalse(maxheight < 2,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many %D subdomains in the vertical directon for mesh height %D", Ndomains, N); 1716af538404SDmitry Karpeev /* Vertical domain limits with an overlap. */ 1717eec7e3faSDmitry Karpeev y[0][0] = PetscMax(ystart - overlap,0); 1718eec7e3faSDmitry Karpeev y[0][1] = PetscMin(ystart + maxheight + overlap,N); 1719eec7e3faSDmitry Karpeev /* Vertical domain limits without an overlap. */ 1720eec7e3faSDmitry Karpeev y[1][0] = ystart; 1721eec7e3faSDmitry Karpeev y[1][1] = PetscMin(ystart + maxheight,N); 1722eec7e3faSDmitry Karpeev xstart = 0; 1723af538404SDmitry Karpeev for (i=0; i<Mdomains; ++i) { 1724af538404SDmitry Karpeev maxwidth = M/Mdomains + ((M % Mdomains) > i); /* Maximal width of subdomain */ 17252c71b3e2SJacob Faibussowitsch PetscCheckFalse(maxwidth < 2,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Too many %D subdomains in the horizontal direction for mesh width %D", Mdomains, M); 1726af538404SDmitry Karpeev /* Horizontal domain limits with an overlap. */ 1727eec7e3faSDmitry Karpeev x[0][0] = PetscMax(xstart - overlap,0); 1728eec7e3faSDmitry Karpeev x[0][1] = PetscMin(xstart + maxwidth + overlap,M); 1729eec7e3faSDmitry Karpeev /* Horizontal domain limits without an overlap. */ 1730eec7e3faSDmitry Karpeev x[1][0] = xstart; 1731eec7e3faSDmitry Karpeev x[1][1] = PetscMin(xstart+maxwidth,M); 17322bbb417fSDmitry Karpeev /* 17332bbb417fSDmitry Karpeev Determine whether this domain intersects this processor's ownership range of pc->pmat. 17342bbb417fSDmitry Karpeev Do this twice: first for the domains with overlaps, and once without. 17352bbb417fSDmitry Karpeev During the first pass create the subcommunicators, and use them on the second pass as well. 17362bbb417fSDmitry Karpeev */ 17372bbb417fSDmitry Karpeev for (q = 0; q < 2; ++q) { 1738cc96b2e5SBarry Smith PetscBool split = PETSC_FALSE; 17392bbb417fSDmitry Karpeev /* 17402bbb417fSDmitry Karpeev domain limits, (xleft, xright) and (ylow, yheigh) are adjusted 17412bbb417fSDmitry Karpeev according to whether the domain with an overlap or without is considered. 17422bbb417fSDmitry Karpeev */ 17432bbb417fSDmitry Karpeev xleft = x[q][0]; xright = x[q][1]; 17442bbb417fSDmitry Karpeev ylow = y[q][0]; yhigh = y[q][1]; 1745c3518bceSDmitry Karpeev PCGASMLocalSubdomainBounds2D(M,N,xleft,ylow,xright,yhigh,first,last,(&xleft_loc),(&ylow_loc),(&xright_loc),(&yhigh_loc),(&nidx)); 1746af538404SDmitry Karpeev nidx *= dof; 1747eec7e3faSDmitry Karpeev n[q] = nidx; 17482bbb417fSDmitry Karpeev /* 1749eec7e3faSDmitry Karpeev Based on the counted number of indices in the local domain *with an overlap*, 1750af538404SDmitry Karpeev construct a subcommunicator of all the processors supporting this domain. 1751eec7e3faSDmitry Karpeev Observe that a domain with an overlap might have nontrivial local support, 1752eec7e3faSDmitry Karpeev while the domain without an overlap might not. Hence, the decision to participate 1753eec7e3faSDmitry Karpeev in the subcommunicator must be based on the domain with an overlap. 17542bbb417fSDmitry Karpeev */ 17552bbb417fSDmitry Karpeev if (q == 0) { 17562fa5cd67SKarl Rupp if (nidx) color = 1; 17572fa5cd67SKarl Rupp else color = MPI_UNDEFINED; 17585f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_split(comm, color, rank, &subcomm)); 1759cc96b2e5SBarry Smith split = PETSC_TRUE; 17602bbb417fSDmitry Karpeev } 1761af538404SDmitry Karpeev /* 1762eec7e3faSDmitry Karpeev Proceed only if the number of local indices *with an overlap* is nonzero. 1763af538404SDmitry Karpeev */ 1764eec7e3faSDmitry Karpeev if (n[0]) { 17652fa5cd67SKarl Rupp if (q == 0) xis = is; 1766af538404SDmitry Karpeev if (q == 1) { 1767af538404SDmitry Karpeev /* 1768eec7e3faSDmitry Karpeev The IS for the no-overlap subdomain shares a communicator with the overlapping domain. 1769af538404SDmitry Karpeev Moreover, if the overlap is zero, the two ISs are identical. 1770af538404SDmitry Karpeev */ 1771b1a0a8a3SJed Brown if (overlap == 0) { 1772eec7e3faSDmitry Karpeev (*is_local)[s] = (*is)[s]; 17735f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectReference((PetscObject)(*is)[s])); 17742bbb417fSDmitry Karpeev continue; 1775d34fcf5fSBarry Smith } else { 17766a4f0f73SDmitry Karpeev xis = is_local; 1777eec7e3faSDmitry Karpeev subcomm = ((PetscObject)(*is)[s])->comm; 17782bbb417fSDmitry Karpeev } 1779af538404SDmitry Karpeev } /* if (q == 1) */ 17800298fd71SBarry Smith idx = NULL; 17815f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(nidx,&idx)); 1782eec7e3faSDmitry Karpeev if (nidx) { 17832bbb417fSDmitry Karpeev k = 0; 17842bbb417fSDmitry Karpeev for (jj=ylow_loc; jj<yhigh_loc; ++jj) { 1785af538404SDmitry Karpeev PetscInt x0 = (jj==ylow_loc) ? xleft_loc : xleft; 1786af538404SDmitry Karpeev PetscInt x1 = (jj==yhigh_loc-1) ? xright_loc : xright; 1787af538404SDmitry Karpeev kk = dof*(M*jj + x0); 17882bbb417fSDmitry Karpeev for (ii=x0; ii<x1; ++ii) { 17892bbb417fSDmitry Karpeev for (d = 0; d < dof; ++d) { 17902bbb417fSDmitry Karpeev idx[k++] = kk++; 1791b1a0a8a3SJed Brown } 1792b1a0a8a3SJed Brown } 1793b1a0a8a3SJed Brown } 1794eec7e3faSDmitry Karpeev } 17955f80ce2aSJacob Faibussowitsch CHKERRQ(ISCreateGeneral(subcomm,nidx,idx,PETSC_OWN_POINTER,(*xis)+s)); 1796cc96b2e5SBarry Smith if (split) { 17975f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_free(&subcomm)); 1798cc96b2e5SBarry Smith } 1799eec7e3faSDmitry Karpeev }/* if (n[0]) */ 18002bbb417fSDmitry Karpeev }/* for (q = 0; q < 2; ++q) */ 18012fa5cd67SKarl Rupp if (n[0]) ++s; 1802af538404SDmitry Karpeev xstart += maxwidth; 1803af538404SDmitry Karpeev } /* for (i = 0; i < Mdomains; ++i) */ 18042bbb417fSDmitry Karpeev ystart += maxheight; 1805af538404SDmitry Karpeev } /* for (j = 0; j < Ndomains; ++j) */ 1806b1a0a8a3SJed Brown PetscFunctionReturn(0); 1807b1a0a8a3SJed Brown } 1808b1a0a8a3SJed Brown 1809b1a0a8a3SJed Brown /*@C 181006b43e7eSDmitry Karpeev PCGASMGetSubdomains - Gets the subdomains supported on this processor 181106b43e7eSDmitry Karpeev for the additive Schwarz preconditioner. 1812b1a0a8a3SJed Brown 1813b1a0a8a3SJed Brown Not Collective 1814b1a0a8a3SJed Brown 1815b1a0a8a3SJed Brown Input Parameter: 1816b1a0a8a3SJed Brown . pc - the preconditioner context 1817b1a0a8a3SJed Brown 1818b1a0a8a3SJed Brown Output Parameters: 1819b1a0a8a3SJed Brown + n - the number of subdomains for this processor (default value = 1) 18200298fd71SBarry Smith . iis - the index sets that define the inner subdomains (without overlap) supported on this processor (can be NULL) 18210298fd71SBarry Smith - ois - the index sets that define the outer subdomains (with overlap) supported on this processor (can be NULL) 1822b1a0a8a3SJed Brown 1823b1a0a8a3SJed Brown Notes: 18246a4f0f73SDmitry Karpeev The user is responsible for destroying the ISs and freeing the returned arrays. 1825b1a0a8a3SJed Brown The IS numbering is in the parallel, global numbering of the vector. 1826b1a0a8a3SJed Brown 1827b1a0a8a3SJed Brown Level: advanced 1828b1a0a8a3SJed Brown 18298f3b4b4dSDmitry Karpeev .seealso: PCGASMSetOverlap(), PCGASMGetSubKSP(), PCGASMCreateSubdomains2D(), 18308f3b4b4dSDmitry Karpeev PCGASMSetSubdomains(), PCGASMGetSubmatrices() 1831b1a0a8a3SJed Brown @*/ 18326a4f0f73SDmitry Karpeev PetscErrorCode PCGASMGetSubdomains(PC pc,PetscInt *n,IS *iis[],IS *ois[]) 1833b1a0a8a3SJed Brown { 1834f746d493SDmitry Karpeev PC_GASM *osm; 1835b1a0a8a3SJed Brown PetscBool match; 18366a4f0f73SDmitry Karpeev PetscInt i; 18375fd66863SKarl Rupp 1838b1a0a8a3SJed Brown PetscFunctionBegin; 1839b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 18405f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match)); 1841*28b400f6SJacob Faibussowitsch PetscCheck(match,PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Incorrect object type: expected %s, got %s instead", PCGASM, ((PetscObject)pc)->type_name); 1842f746d493SDmitry Karpeev osm = (PC_GASM*)pc->data; 1843ab3e923fSDmitry Karpeev if (n) *n = osm->n; 18446a4f0f73SDmitry Karpeev if (iis) { 18455f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(osm->n, iis)); 18466a4f0f73SDmitry Karpeev } 18476a4f0f73SDmitry Karpeev if (ois) { 18485f80ce2aSJacob Faibussowitsch CHKERRQ(PetscMalloc1(osm->n, ois)); 18496a4f0f73SDmitry Karpeev } 18506a4f0f73SDmitry Karpeev if (iis || ois) { 18516a4f0f73SDmitry Karpeev for (i = 0; i < osm->n; ++i) { 18526a4f0f73SDmitry Karpeev if (iis) (*iis)[i] = osm->iis[i]; 18536a4f0f73SDmitry Karpeev if (ois) (*ois)[i] = osm->ois[i]; 18546a4f0f73SDmitry Karpeev } 1855b1a0a8a3SJed Brown } 1856b1a0a8a3SJed Brown PetscFunctionReturn(0); 1857b1a0a8a3SJed Brown } 1858b1a0a8a3SJed Brown 1859b1a0a8a3SJed Brown /*@C 186006b43e7eSDmitry Karpeev PCGASMGetSubmatrices - Gets the local submatrices (for this processor 1861b1a0a8a3SJed Brown only) for the additive Schwarz preconditioner. 1862b1a0a8a3SJed Brown 1863b1a0a8a3SJed Brown Not Collective 1864b1a0a8a3SJed Brown 1865b1a0a8a3SJed Brown Input Parameter: 1866b1a0a8a3SJed Brown . pc - the preconditioner context 1867b1a0a8a3SJed Brown 1868b1a0a8a3SJed Brown Output Parameters: 1869b1a0a8a3SJed Brown + n - the number of matrices for this processor (default value = 1) 1870b1a0a8a3SJed Brown - mat - the matrices 1871b1a0a8a3SJed Brown 187295452b02SPatrick Sanan Notes: 187395452b02SPatrick Sanan matrices returned by this routine have the same communicators as the index sets (IS) 18748f3b4b4dSDmitry Karpeev used to define subdomains in PCGASMSetSubdomains() 1875b1a0a8a3SJed Brown Level: advanced 1876b1a0a8a3SJed Brown 18778f3b4b4dSDmitry Karpeev .seealso: PCGASMSetOverlap(), PCGASMGetSubKSP(), 187806b43e7eSDmitry Karpeev PCGASMCreateSubdomains2D(), PCGASMSetSubdomains(), PCGASMGetSubdomains() 1879b1a0a8a3SJed Brown @*/ 188006b43e7eSDmitry Karpeev PetscErrorCode PCGASMGetSubmatrices(PC pc,PetscInt *n,Mat *mat[]) 1881b1a0a8a3SJed Brown { 1882f746d493SDmitry Karpeev PC_GASM *osm; 1883b1a0a8a3SJed Brown PetscBool match; 1884b1a0a8a3SJed Brown 1885b1a0a8a3SJed Brown PetscFunctionBegin; 1886b1a0a8a3SJed Brown PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1887b1a0a8a3SJed Brown PetscValidIntPointer(n,2); 1888b1a0a8a3SJed Brown if (mat) PetscValidPointer(mat,3); 1889*28b400f6SJacob Faibussowitsch PetscCheck(pc->setupcalled,PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONGSTATE,"Must call after KSPSetUp() or PCSetUp()."); 18905f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match)); 1891*28b400f6SJacob Faibussowitsch PetscCheck(match,PetscObjectComm((PetscObject)pc), PETSC_ERR_ARG_WRONG, "Expected %s, got %s instead", PCGASM, ((PetscObject)pc)->type_name); 1892f746d493SDmitry Karpeev osm = (PC_GASM*)pc->data; 1893ab3e923fSDmitry Karpeev if (n) *n = osm->n; 1894b1a0a8a3SJed Brown if (mat) *mat = osm->pmat; 1895b1a0a8a3SJed Brown PetscFunctionReturn(0); 1896b1a0a8a3SJed Brown } 1897d709ea83SDmitry Karpeev 1898d709ea83SDmitry Karpeev /*@ 18998f3b4b4dSDmitry Karpeev PCGASMSetUseDMSubdomains - Indicates whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible. 1900d709ea83SDmitry Karpeev Logically Collective 1901d709ea83SDmitry Karpeev 1902d8d19677SJose E. Roman Input Parameters: 1903d709ea83SDmitry Karpeev + pc - the preconditioner 1904d709ea83SDmitry Karpeev - flg - boolean indicating whether to use subdomains defined by the DM 1905d709ea83SDmitry Karpeev 1906d709ea83SDmitry Karpeev Options Database Key: 19078f3b4b4dSDmitry Karpeev . -pc_gasm_dm_subdomains -pc_gasm_overlap -pc_gasm_total_subdomains 1908d709ea83SDmitry Karpeev 1909d709ea83SDmitry Karpeev Level: intermediate 1910d709ea83SDmitry Karpeev 1911d709ea83SDmitry Karpeev Notes: 19128f3b4b4dSDmitry Karpeev PCGASMSetSubdomains(), PCGASMSetTotalSubdomains() or PCGASMSetOverlap() take precedence over PCGASMSetUseDMSubdomains(), 19138f3b4b4dSDmitry Karpeev so setting PCGASMSetSubdomains() with nontrivial subdomain ISs or any of PCGASMSetTotalSubdomains() and PCGASMSetOverlap() 19148f3b4b4dSDmitry Karpeev automatically turns the latter off. 1915d709ea83SDmitry Karpeev 19168f3b4b4dSDmitry Karpeev .seealso: PCGASMGetUseDMSubdomains(), PCGASMSetSubdomains(), PCGASMSetOverlap() 1917d709ea83SDmitry Karpeev PCGASMCreateSubdomains2D() 1918d709ea83SDmitry Karpeev @*/ 19198f3b4b4dSDmitry Karpeev PetscErrorCode PCGASMSetUseDMSubdomains(PC pc,PetscBool flg) 1920d709ea83SDmitry Karpeev { 1921d709ea83SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 1922d709ea83SDmitry Karpeev PetscBool match; 1923d709ea83SDmitry Karpeev 1924d709ea83SDmitry Karpeev PetscFunctionBegin; 1925d709ea83SDmitry Karpeev PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1926d709ea83SDmitry Karpeev PetscValidLogicalCollectiveBool(pc,flg,2); 1927*28b400f6SJacob Faibussowitsch PetscCheck(!pc->setupcalled,((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONGSTATE,"Not for a setup PC."); 19285f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match)); 1929d709ea83SDmitry Karpeev if (match) { 19308f3b4b4dSDmitry Karpeev if (!osm->user_subdomains && osm->N == PETSC_DETERMINE && osm->overlap < 0) { 1931d709ea83SDmitry Karpeev osm->dm_subdomains = flg; 1932d709ea83SDmitry Karpeev } 19338f3b4b4dSDmitry Karpeev } 1934d709ea83SDmitry Karpeev PetscFunctionReturn(0); 1935d709ea83SDmitry Karpeev } 1936d709ea83SDmitry Karpeev 1937d709ea83SDmitry Karpeev /*@ 19388f3b4b4dSDmitry Karpeev PCGASMGetUseDMSubdomains - Returns flag indicating whether to use DMCreateDomainDecomposition() to define the subdomains, whenever possible. 1939d709ea83SDmitry Karpeev Not Collective 1940d709ea83SDmitry Karpeev 1941d709ea83SDmitry Karpeev Input Parameter: 1942d709ea83SDmitry Karpeev . pc - the preconditioner 1943d709ea83SDmitry Karpeev 1944d709ea83SDmitry Karpeev Output Parameter: 1945d709ea83SDmitry Karpeev . flg - boolean indicating whether to use subdomains defined by the DM 1946d709ea83SDmitry Karpeev 1947d709ea83SDmitry Karpeev Level: intermediate 1948d709ea83SDmitry Karpeev 19498f3b4b4dSDmitry Karpeev .seealso: PCGASMSetUseDMSubdomains(), PCGASMSetOverlap() 1950d709ea83SDmitry Karpeev PCGASMCreateSubdomains2D() 1951d709ea83SDmitry Karpeev @*/ 19528f3b4b4dSDmitry Karpeev PetscErrorCode PCGASMGetUseDMSubdomains(PC pc,PetscBool* flg) 1953d709ea83SDmitry Karpeev { 1954d709ea83SDmitry Karpeev PC_GASM *osm = (PC_GASM*)pc->data; 1955d709ea83SDmitry Karpeev PetscBool match; 1956d709ea83SDmitry Karpeev 1957d709ea83SDmitry Karpeev PetscFunctionBegin; 1958d709ea83SDmitry Karpeev PetscValidHeaderSpecific(pc,PC_CLASSID,1); 1959534a8f05SLisandro Dalcin PetscValidBoolPointer(flg,2); 19605f80ce2aSJacob Faibussowitsch CHKERRQ(PetscObjectTypeCompare((PetscObject)pc,PCGASM,&match)); 1961d709ea83SDmitry Karpeev if (match) { 1962d709ea83SDmitry Karpeev if (flg) *flg = osm->dm_subdomains; 1963d709ea83SDmitry Karpeev } 1964d709ea83SDmitry Karpeev PetscFunctionReturn(0); 1965d709ea83SDmitry Karpeev } 1966