1dba47a55SKris Buschelman 24b9ad928SBarry Smith /* 33f457be1SHong Zhang This file defines a "solve the problem redundantly on each subgroup of processor" preconditioner. 44b9ad928SBarry Smith */ 5c6db04a5SJed Brown #include <private/pcimpl.h> /*I "petscpc.h" I*/ 6c6db04a5SJed Brown #include <petscksp.h> /*I "petscksp.h" I*/ 74b9ad928SBarry Smith 84b9ad928SBarry Smith typedef struct { 93e065800SHong Zhang KSP ksp; 104b9ad928SBarry Smith PC pc; /* actual preconditioner used on each processor */ 117adad957SLisandro Dalcin Vec xsub,ysub; /* vectors of a subcommunicator to hold parallel vectors of ((PetscObject)pc)->comm */ 123f457be1SHong Zhang Vec xdup,ydup; /* parallel vector that congregates xsub or ysub facilitating vector scattering */ 13b3804887SHong Zhang Mat pmats; /* matrix and optional preconditioner matrix belong to a subcommunicator */ 143f457be1SHong Zhang VecScatter scatterin,scatterout; /* scatter used to move all values to each processor group (subcommunicator) */ 15ace3abfcSBarry Smith PetscBool useparallelmat; 16c540e29cSHong Zhang PetscSubcomm psubcomm; 171fbd8f88SHong Zhang PetscInt nsubcomm; /* num of data structure PetscSubcomm */ 184b9ad928SBarry Smith } PC_Redundant; 194b9ad928SBarry Smith 204b9ad928SBarry Smith #undef __FUNCT__ 214b9ad928SBarry Smith #define __FUNCT__ "PCView_Redundant" 226849ba73SBarry Smith static PetscErrorCode PCView_Redundant(PC pc,PetscViewer viewer) 234b9ad928SBarry Smith { 244b9ad928SBarry Smith PC_Redundant *red = (PC_Redundant*)pc->data; 25dfbe8321SBarry Smith PetscErrorCode ierr; 26ace3abfcSBarry Smith PetscBool iascii,isstring; 2703ccd0b4SBarry Smith PetscViewer subviewer; 284b9ad928SBarry Smith 294b9ad928SBarry Smith PetscFunctionBegin; 302692d6eeSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 312692d6eeSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr); 3232077d6dSBarry Smith if (iascii) { 3303ccd0b4SBarry Smith if (!red->psubcomm) { 3403ccd0b4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Redundant preconditioner: Not yet setup\n");CHKERRQ(ierr); 3503ccd0b4SBarry Smith } else { 363e065800SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Redundant preconditioner: First (color=0) of %D PCs follows\n",red->nsubcomm);CHKERRQ(ierr); 377adad957SLisandro Dalcin ierr = PetscViewerGetSubcomm(viewer,((PetscObject)red->pc)->comm,&subviewer);CHKERRQ(ierr); 38f5dd71faSBarry Smith if (!red->psubcomm->color) { /* only view first redundant pc */ 394b9ad928SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 403e065800SHong Zhang ierr = KSPView(red->ksp,subviewer);CHKERRQ(ierr); 414b9ad928SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 424b9ad928SBarry Smith } 437adad957SLisandro Dalcin ierr = PetscViewerRestoreSubcomm(viewer,((PetscObject)red->pc)->comm,&subviewer);CHKERRQ(ierr); 444b9ad928SBarry Smith } 4503ccd0b4SBarry Smith } else if (isstring) { 4603ccd0b4SBarry Smith ierr = PetscViewerStringSPrintf(viewer," Redundant solver preconditioner");CHKERRQ(ierr); 474b9ad928SBarry Smith } else { 4865e19b50SBarry Smith SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PC redundant",((PetscObject)viewer)->type_name); 494b9ad928SBarry Smith } 504b9ad928SBarry Smith PetscFunctionReturn(0); 514b9ad928SBarry Smith } 524b9ad928SBarry Smith 534b9ad928SBarry Smith #undef __FUNCT__ 544b9ad928SBarry Smith #define __FUNCT__ "PCSetUp_Redundant" 556849ba73SBarry Smith static PetscErrorCode PCSetUp_Redundant(PC pc) 564b9ad928SBarry Smith { 574b9ad928SBarry Smith PC_Redundant *red = (PC_Redundant*)pc->data; 58dfbe8321SBarry Smith PetscErrorCode ierr; 5945fc02eaSBarry Smith PetscInt mstart,mend,mlocal,m,mlocal_sub,rstart_sub,rend_sub,mloc_sub; 6013f74950SBarry Smith PetscMPIInt size; 614b9ad928SBarry Smith MatReuse reuse = MAT_INITIAL_MATRIX; 624b9ad928SBarry Smith MatStructure str = DIFFERENT_NONZERO_PATTERN; 637adad957SLisandro Dalcin MPI_Comm comm = ((PetscObject)pc)->comm,subcomm; 6423ce1328SBarry Smith Vec vec; 653f457be1SHong Zhang PetscMPIInt subsize,subrank; 661fbd8f88SHong Zhang const char *prefix; 67b862ddfaSBarry Smith const PetscInt *range; 683f457be1SHong Zhang 694b9ad928SBarry Smith PetscFunctionBegin; 7023ce1328SBarry Smith ierr = MatGetVecs(pc->pmat,&vec,0);CHKERRQ(ierr); 7123ce1328SBarry Smith ierr = VecGetSize(vec,&m);CHKERRQ(ierr); 721fbd8f88SHong Zhang 734b9ad928SBarry Smith if (!pc->setupcalled) { 745f06b7aaSBarry Smith if (!red->psubcomm) { 75d8a68f86SHong Zhang ierr = PetscSubcommCreate(comm,&red->psubcomm);CHKERRQ(ierr); 76d8a68f86SHong Zhang ierr = PetscSubcommSetNumber(red->psubcomm,red->nsubcomm);CHKERRQ(ierr); 77d8a68f86SHong Zhang ierr = PetscSubcommSetType(red->psubcomm,PETSC_SUBCOMM_INTERLACED);CHKERRQ(ierr); 781fbd8f88SHong Zhang ierr = PetscLogObjectMemory(pc,sizeof(PetscSubcomm));CHKERRQ(ierr); 791fbd8f88SHong Zhang 801fbd8f88SHong Zhang /* create a new PC that processors in each subcomm have copy of */ 810d7810c8SBarry Smith subcomm = red->psubcomm->comm; 825f06b7aaSBarry Smith ierr = KSPCreate(subcomm,&red->ksp);CHKERRQ(ierr); 835f06b7aaSBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)red->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 845f06b7aaSBarry Smith ierr = PetscLogObjectParent(pc,red->ksp);CHKERRQ(ierr); 855f06b7aaSBarry Smith ierr = KSPSetType(red->ksp,KSPPREONLY);CHKERRQ(ierr); 865f06b7aaSBarry Smith ierr = KSPGetPC(red->ksp,&red->pc);CHKERRQ(ierr); 87cf52b8b1SHong Zhang ierr = PCSetType(red->pc,PCLU);CHKERRQ(ierr); 88cf52b8b1SHong Zhang 891fbd8f88SHong Zhang ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 905f06b7aaSBarry Smith ierr = KSPSetOptionsPrefix(red->ksp,prefix);CHKERRQ(ierr); 915f06b7aaSBarry Smith ierr = KSPAppendOptionsPrefix(red->ksp,"redundant_");CHKERRQ(ierr); 925f06b7aaSBarry Smith } else { 935f06b7aaSBarry Smith subcomm = red->psubcomm->comm; 945f06b7aaSBarry Smith } 951fbd8f88SHong Zhang 963f457be1SHong Zhang /* create working vectors xsub/ysub and xdup/ydup */ 9723ce1328SBarry Smith ierr = VecGetLocalSize(vec,&mlocal);CHKERRQ(ierr); 983f457be1SHong Zhang ierr = VecGetOwnershipRange(vec,&mstart,&mend);CHKERRQ(ierr); 994b9ad928SBarry Smith 1003f457be1SHong Zhang /* get local size of xsub/ysub */ 1011fbd8f88SHong Zhang ierr = MPI_Comm_size(subcomm,&subsize);CHKERRQ(ierr); 1021fbd8f88SHong Zhang ierr = MPI_Comm_rank(subcomm,&subrank);CHKERRQ(ierr); 103b862ddfaSBarry Smith ierr = MatGetOwnershipRanges(pc->pmat,&range);CHKERRQ(ierr); 104b862ddfaSBarry Smith rstart_sub = range[red->psubcomm->n*subrank]; /* rstart in xsub/ysub */ 1053f457be1SHong Zhang if (subrank+1 < subsize){ 106b862ddfaSBarry Smith rend_sub = range[red->psubcomm->n*(subrank+1)]; 1073f457be1SHong Zhang } else { 1083f457be1SHong Zhang rend_sub = m; 1093f457be1SHong Zhang } 1103f457be1SHong Zhang mloc_sub = rend_sub - rstart_sub; 1111fbd8f88SHong Zhang ierr = VecCreateMPI(subcomm,mloc_sub,PETSC_DECIDE,&red->ysub);CHKERRQ(ierr); 1123f457be1SHong Zhang /* create xsub with empty local arrays, because xdup's arrays will be placed into it */ 1131fbd8f88SHong Zhang ierr = VecCreateMPIWithArray(subcomm,mloc_sub,PETSC_DECIDE,PETSC_NULL,&red->xsub);CHKERRQ(ierr); 1143f457be1SHong Zhang 1153f457be1SHong Zhang /* create xdup and ydup. ydup has empty local arrays because ysub's arrays will be place into it. 1167adad957SLisandro Dalcin Note: we use communicator dupcomm, not ((PetscObject)pc)->comm! */ 1171fbd8f88SHong Zhang ierr = VecCreateMPI(red->psubcomm->dupparent,mloc_sub,PETSC_DECIDE,&red->xdup);CHKERRQ(ierr); 1181fbd8f88SHong Zhang ierr = VecCreateMPIWithArray(red->psubcomm->dupparent,mloc_sub,PETSC_DECIDE,PETSC_NULL,&red->ydup);CHKERRQ(ierr); 1193f457be1SHong Zhang 1203f457be1SHong Zhang /* create vec scatters */ 1213f457be1SHong Zhang if (!red->scatterin){ 1223f457be1SHong Zhang IS is1,is2; 1233f457be1SHong Zhang PetscInt *idx1,*idx2,i,j,k; 12445fc02eaSBarry Smith 1251d79065fSBarry Smith ierr = PetscMalloc2(red->psubcomm->n*mlocal,PetscInt,&idx1,red->psubcomm->n*mlocal,PetscInt,&idx2);CHKERRQ(ierr); 1263f457be1SHong Zhang j = 0; 1271fbd8f88SHong Zhang for (k=0; k<red->psubcomm->n; k++){ 1283f457be1SHong Zhang for (i=mstart; i<mend; i++){ 1293f457be1SHong Zhang idx1[j] = i; 1303f457be1SHong Zhang idx2[j++] = i + m*k; 1313f457be1SHong Zhang } 1323f457be1SHong Zhang } 13370b3c8c7SBarry Smith ierr = ISCreateGeneral(comm,red->psubcomm->n*mlocal,idx1,PETSC_COPY_VALUES,&is1);CHKERRQ(ierr); 13470b3c8c7SBarry Smith ierr = ISCreateGeneral(comm,red->psubcomm->n*mlocal,idx2,PETSC_COPY_VALUES,&is2);CHKERRQ(ierr); 1353f457be1SHong Zhang ierr = VecScatterCreate(vec,is1,red->xdup,is2,&red->scatterin);CHKERRQ(ierr); 136fcfd50ebSBarry Smith ierr = ISDestroy(&is1);CHKERRQ(ierr); 137fcfd50ebSBarry Smith ierr = ISDestroy(&is2);CHKERRQ(ierr); 1383f457be1SHong Zhang 1391fbd8f88SHong Zhang ierr = ISCreateStride(comm,mlocal,mstart+ red->psubcomm->color*m,1,&is1);CHKERRQ(ierr); 1403f457be1SHong Zhang ierr = ISCreateStride(comm,mlocal,mstart,1,&is2);CHKERRQ(ierr); 1413f457be1SHong Zhang ierr = VecScatterCreate(red->xdup,is1,vec,is2,&red->scatterout);CHKERRQ(ierr); 142fcfd50ebSBarry Smith ierr = ISDestroy(&is1);CHKERRQ(ierr); 143fcfd50ebSBarry Smith ierr = ISDestroy(&is2);CHKERRQ(ierr); 1441d79065fSBarry Smith ierr = PetscFree2(idx1,idx2);CHKERRQ(ierr); 1454b9ad928SBarry Smith } 1464b9ad928SBarry Smith } 1476bf464f9SBarry Smith ierr = VecDestroy(&vec);CHKERRQ(ierr); 1484b9ad928SBarry Smith 1494b9ad928SBarry Smith /* if pmatrix set by user is sequential then we do not need to gather the parallel matrix */ 1503f457be1SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1514b9ad928SBarry Smith if (size == 1) { 1524b9ad928SBarry Smith red->useparallelmat = PETSC_FALSE; 1534b9ad928SBarry Smith } 1544b9ad928SBarry Smith 1554b9ad928SBarry Smith if (red->useparallelmat) { 1564b9ad928SBarry Smith if (pc->setupcalled == 1 && pc->flag == DIFFERENT_NONZERO_PATTERN) { 1574b9ad928SBarry Smith /* destroy old matrices */ 1586bf464f9SBarry Smith ierr = MatDestroy(&red->pmats);CHKERRQ(ierr); 1594b9ad928SBarry Smith } else if (pc->setupcalled == 1) { 1604b9ad928SBarry Smith reuse = MAT_REUSE_MATRIX; 1614b9ad928SBarry Smith str = SAME_NONZERO_PATTERN; 1624b9ad928SBarry Smith } 1634b9ad928SBarry Smith 1643f457be1SHong Zhang /* grab the parallel matrix and put it into processors of a subcomminicator */ 165f664ae05SHong Zhang /*--------------------------------------------------------------------------*/ 166f664ae05SHong Zhang ierr = VecGetLocalSize(red->ysub,&mlocal_sub);CHKERRQ(ierr); 16769db28dcSHong Zhang ierr = MatGetRedundantMatrix(pc->pmat,red->psubcomm->n,red->psubcomm->comm,mlocal_sub,reuse,&red->pmats);CHKERRQ(ierr); 1683f457be1SHong Zhang /* tell PC of the subcommunicator its operators */ 16990f1c854SHong Zhang ierr = KSPSetOperators(red->ksp,red->pmats,red->pmats,str);CHKERRQ(ierr); 1704b9ad928SBarry Smith } else { 17190f1c854SHong Zhang ierr = KSPSetOperators(red->ksp,pc->mat,pc->pmat,pc->flag);CHKERRQ(ierr); 1724b9ad928SBarry Smith } 1730c24e6a1SHong Zhang if (pc->setfromoptionscalled){ 1743e065800SHong Zhang ierr = KSPSetFromOptions(red->ksp);CHKERRQ(ierr); 1750c24e6a1SHong Zhang } 1763e065800SHong Zhang ierr = KSPSetUp(red->ksp);CHKERRQ(ierr); 1774b9ad928SBarry Smith PetscFunctionReturn(0); 1784b9ad928SBarry Smith } 1794b9ad928SBarry Smith 1804b9ad928SBarry Smith #undef __FUNCT__ 1814b9ad928SBarry Smith #define __FUNCT__ "PCApply_Redundant" 1826849ba73SBarry Smith static PetscErrorCode PCApply_Redundant(PC pc,Vec x,Vec y) 1834b9ad928SBarry Smith { 1844b9ad928SBarry Smith PC_Redundant *red = (PC_Redundant*)pc->data; 185dfbe8321SBarry Smith PetscErrorCode ierr; 1863f457be1SHong Zhang PetscScalar *array; 1874b9ad928SBarry Smith 1884b9ad928SBarry Smith PetscFunctionBegin; 1893f457be1SHong Zhang /* scatter x to xdup */ 190ca9f406cSSatish Balay ierr = VecScatterBegin(red->scatterin,x,red->xdup,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 191ca9f406cSSatish Balay ierr = VecScatterEnd(red->scatterin,x,red->xdup,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1923f457be1SHong Zhang 1933f457be1SHong Zhang /* place xdup's local array into xsub */ 1943f457be1SHong Zhang ierr = VecGetArray(red->xdup,&array);CHKERRQ(ierr); 1953f457be1SHong Zhang ierr = VecPlaceArray(red->xsub,(const PetscScalar*)array);CHKERRQ(ierr); 1964b9ad928SBarry Smith 1974b9ad928SBarry Smith /* apply preconditioner on each processor */ 19883ab6a24SBarry Smith ierr = KSPSolve(red->ksp,red->xsub,red->ysub);CHKERRQ(ierr); 1993f457be1SHong Zhang ierr = VecResetArray(red->xsub);CHKERRQ(ierr); 2003f457be1SHong Zhang ierr = VecRestoreArray(red->xdup,&array);CHKERRQ(ierr); 2014b9ad928SBarry Smith 2023f457be1SHong Zhang /* place ysub's local array into ydup */ 2033f457be1SHong Zhang ierr = VecGetArray(red->ysub,&array);CHKERRQ(ierr); 2043f457be1SHong Zhang ierr = VecPlaceArray(red->ydup,(const PetscScalar*)array);CHKERRQ(ierr); 2053f457be1SHong Zhang 2063f457be1SHong Zhang /* scatter ydup to y */ 207ca9f406cSSatish Balay ierr = VecScatterBegin(red->scatterout,red->ydup,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 208ca9f406cSSatish Balay ierr = VecScatterEnd(red->scatterout,red->ydup,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 2093f457be1SHong Zhang ierr = VecResetArray(red->ydup);CHKERRQ(ierr); 2103f457be1SHong Zhang ierr = VecRestoreArray(red->ysub,&array);CHKERRQ(ierr); 2114b9ad928SBarry Smith PetscFunctionReturn(0); 2124b9ad928SBarry Smith } 2134b9ad928SBarry Smith 2144b9ad928SBarry Smith #undef __FUNCT__ 2151ea5a559SBarry Smith #define __FUNCT__ "PCReset_Redundant" 2161ea5a559SBarry Smith static PetscErrorCode PCReset_Redundant(PC pc) 2174b9ad928SBarry Smith { 2184b9ad928SBarry Smith PC_Redundant *red = (PC_Redundant*)pc->data; 219dfbe8321SBarry Smith PetscErrorCode ierr; 2204b9ad928SBarry Smith 2214b9ad928SBarry Smith PetscFunctionBegin; 2226bf464f9SBarry Smith ierr = VecScatterDestroy(&red->scatterin);CHKERRQ(ierr); 2236bf464f9SBarry Smith ierr = VecScatterDestroy(&red->scatterout);CHKERRQ(ierr); 2246bf464f9SBarry Smith ierr = VecDestroy(&red->ysub);CHKERRQ(ierr); 2256bf464f9SBarry Smith ierr = VecDestroy(&red->xsub);CHKERRQ(ierr); 2266bf464f9SBarry Smith ierr = VecDestroy(&red->xdup);CHKERRQ(ierr); 2276bf464f9SBarry Smith ierr = VecDestroy(&red->ydup);CHKERRQ(ierr); 2286bf464f9SBarry Smith ierr = MatDestroy(&red->pmats);CHKERRQ(ierr); 229*546c30cbSLisandro Dalcin if (red->ksp) {ierr = KSPReset(red->ksp);CHKERRQ(ierr);} 2301ea5a559SBarry Smith PetscFunctionReturn(0); 2311ea5a559SBarry Smith } 2321ea5a559SBarry Smith 2331ea5a559SBarry Smith #undef __FUNCT__ 2341ea5a559SBarry Smith #define __FUNCT__ "PCDestroy_Redundant" 2351ea5a559SBarry Smith static PetscErrorCode PCDestroy_Redundant(PC pc) 2361ea5a559SBarry Smith { 2371ea5a559SBarry Smith PC_Redundant *red = (PC_Redundant*)pc->data; 2381ea5a559SBarry Smith PetscErrorCode ierr; 2391ea5a559SBarry Smith 2401ea5a559SBarry Smith PetscFunctionBegin; 2411ea5a559SBarry Smith ierr = PCReset_Redundant(pc);CHKERRQ(ierr); 2426bf464f9SBarry Smith ierr = KSPDestroy(&red->ksp);CHKERRQ(ierr); 2436bf464f9SBarry Smith ierr = PetscSubcommDestroy(&red->psubcomm);CHKERRQ(ierr); 244c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 2454b9ad928SBarry Smith PetscFunctionReturn(0); 2464b9ad928SBarry Smith } 2474b9ad928SBarry Smith 2484b9ad928SBarry Smith #undef __FUNCT__ 2494b9ad928SBarry Smith #define __FUNCT__ "PCSetFromOptions_Redundant" 2506849ba73SBarry Smith static PetscErrorCode PCSetFromOptions_Redundant(PC pc) 2514b9ad928SBarry Smith { 252a98ce0f4SHong Zhang PetscErrorCode ierr; 253a98ce0f4SHong Zhang PC_Redundant *red = (PC_Redundant*)pc->data; 254a98ce0f4SHong Zhang 2554b9ad928SBarry Smith PetscFunctionBegin; 256a98ce0f4SHong Zhang ierr = PetscOptionsHead("Redundant options");CHKERRQ(ierr); 25709a6bc64SHong Zhang ierr = PetscOptionsInt("-pc_redundant_number","Number of redundant pc","PCRedundantSetNumber",red->nsubcomm,&red->nsubcomm,0);CHKERRQ(ierr); 258a98ce0f4SHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 2594b9ad928SBarry Smith PetscFunctionReturn(0); 2604b9ad928SBarry Smith } 2614b9ad928SBarry Smith 2624b9ad928SBarry Smith EXTERN_C_BEGIN 2634b9ad928SBarry Smith #undef __FUNCT__ 26409a6bc64SHong Zhang #define __FUNCT__ "PCRedundantSetNumber_Redundant" 2657087cfbeSBarry Smith PetscErrorCode PCRedundantSetNumber_Redundant(PC pc,PetscInt nreds) 26609a6bc64SHong Zhang { 26709a6bc64SHong Zhang PC_Redundant *red = (PC_Redundant*)pc->data; 26809a6bc64SHong Zhang 26909a6bc64SHong Zhang PetscFunctionBegin; 27009a6bc64SHong Zhang red->nsubcomm = nreds; 27109a6bc64SHong Zhang PetscFunctionReturn(0); 27209a6bc64SHong Zhang } 27309a6bc64SHong Zhang EXTERN_C_END 27409a6bc64SHong Zhang 27509a6bc64SHong Zhang #undef __FUNCT__ 27609a6bc64SHong Zhang #define __FUNCT__ "PCRedundantSetNumber" 27709a6bc64SHong Zhang /*@ 27809a6bc64SHong Zhang PCRedundantSetNumber - Sets the number of redundant preconditioner contexts. 27909a6bc64SHong Zhang 2803f9fe445SBarry Smith Logically Collective on PC 28109a6bc64SHong Zhang 28209a6bc64SHong Zhang Input Parameters: 28309a6bc64SHong Zhang + pc - the preconditioner context 2849b21d695SBarry Smith - nredundant - number of redundant preconditioner contexts; for example if you are using 64 MPI processes and 2859b21d695SBarry Smith use an nredundant of 4 there will be 4 parallel solves each on 16 = 64/4 processes. 28609a6bc64SHong Zhang 28709a6bc64SHong Zhang Level: advanced 28809a6bc64SHong Zhang 28909a6bc64SHong Zhang .keywords: PC, redundant solve 29009a6bc64SHong Zhang @*/ 2917087cfbeSBarry Smith PetscErrorCode PCRedundantSetNumber(PC pc,PetscInt nredundant) 29209a6bc64SHong Zhang { 2934ac538c5SBarry Smith PetscErrorCode ierr; 29409a6bc64SHong Zhang 29509a6bc64SHong Zhang PetscFunctionBegin; 2960700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 29765e19b50SBarry Smith if (nredundant <= 0) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG, "num of redundant pc %D must be positive",nredundant); 2984ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCRedundantSetNumber_C",(PC,PetscInt),(pc,nredundant));CHKERRQ(ierr); 29909a6bc64SHong Zhang PetscFunctionReturn(0); 30009a6bc64SHong Zhang } 30109a6bc64SHong Zhang 30209a6bc64SHong Zhang EXTERN_C_BEGIN 30309a6bc64SHong Zhang #undef __FUNCT__ 3044b9ad928SBarry Smith #define __FUNCT__ "PCRedundantSetScatter_Redundant" 3057087cfbeSBarry Smith PetscErrorCode PCRedundantSetScatter_Redundant(PC pc,VecScatter in,VecScatter out) 3064b9ad928SBarry Smith { 3074b9ad928SBarry Smith PC_Redundant *red = (PC_Redundant*)pc->data; 308dfbe8321SBarry Smith PetscErrorCode ierr; 3094b9ad928SBarry Smith 3104b9ad928SBarry Smith PetscFunctionBegin; 3114b9ad928SBarry Smith ierr = PetscObjectReference((PetscObject)in);CHKERRQ(ierr); 3126bf464f9SBarry Smith ierr = VecScatterDestroy(&red->scatterin);CHKERRQ(ierr); 313c3122656SLisandro Dalcin red->scatterin = in; 3144b9ad928SBarry Smith ierr = PetscObjectReference((PetscObject)out);CHKERRQ(ierr); 3156bf464f9SBarry Smith ierr = VecScatterDestroy(&red->scatterout);CHKERRQ(ierr); 316c3122656SLisandro Dalcin red->scatterout = out; 3174b9ad928SBarry Smith PetscFunctionReturn(0); 3184b9ad928SBarry Smith } 3194b9ad928SBarry Smith EXTERN_C_END 3204b9ad928SBarry Smith 3214b9ad928SBarry Smith #undef __FUNCT__ 3224b9ad928SBarry Smith #define __FUNCT__ "PCRedundantSetScatter" 3234b9ad928SBarry Smith /*@ 3244b9ad928SBarry Smith PCRedundantSetScatter - Sets the scatter used to copy values into the 3254b9ad928SBarry Smith redundant local solve and the scatter to move them back into the global 3264b9ad928SBarry Smith vector. 3274b9ad928SBarry Smith 3283f9fe445SBarry Smith Logically Collective on PC 3294b9ad928SBarry Smith 3304b9ad928SBarry Smith Input Parameters: 3314b9ad928SBarry Smith + pc - the preconditioner context 3324b9ad928SBarry Smith . in - the scatter to move the values in 3334b9ad928SBarry Smith - out - the scatter to move them out 3344b9ad928SBarry Smith 3354b9ad928SBarry Smith Level: advanced 3364b9ad928SBarry Smith 3374b9ad928SBarry Smith .keywords: PC, redundant solve 3384b9ad928SBarry Smith @*/ 3397087cfbeSBarry Smith PetscErrorCode PCRedundantSetScatter(PC pc,VecScatter in,VecScatter out) 3404b9ad928SBarry Smith { 3414ac538c5SBarry Smith PetscErrorCode ierr; 3424b9ad928SBarry Smith 3434b9ad928SBarry Smith PetscFunctionBegin; 3440700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 3450700a824SBarry Smith PetscValidHeaderSpecific(in,VEC_SCATTER_CLASSID,2); 3460700a824SBarry Smith PetscValidHeaderSpecific(out,VEC_SCATTER_CLASSID,3); 3474ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCRedundantSetScatter_C",(PC,VecScatter,VecScatter),(pc,in,out));CHKERRQ(ierr); 3484b9ad928SBarry Smith PetscFunctionReturn(0); 3494b9ad928SBarry Smith } 3504b9ad928SBarry Smith 3514b9ad928SBarry Smith EXTERN_C_BEGIN 3524b9ad928SBarry Smith #undef __FUNCT__ 35383ab6a24SBarry Smith #define __FUNCT__ "PCRedundantGetKSP_Redundant" 35483ab6a24SBarry Smith PetscErrorCode PCRedundantGetKSP_Redundant(PC pc,KSP *innerksp) 3554b9ad928SBarry Smith { 3565f06b7aaSBarry Smith PetscErrorCode ierr; 3574b9ad928SBarry Smith PC_Redundant *red = (PC_Redundant*)pc->data; 3585f06b7aaSBarry Smith MPI_Comm comm,subcomm; 3595f06b7aaSBarry Smith const char *prefix; 3604b9ad928SBarry Smith 3614b9ad928SBarry Smith PetscFunctionBegin; 3625f06b7aaSBarry Smith if (!red->psubcomm) { 3635f06b7aaSBarry Smith ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 364d8a68f86SHong Zhang ierr = PetscSubcommCreate(comm,&red->psubcomm);CHKERRQ(ierr); 365d8a68f86SHong Zhang ierr = PetscSubcommSetNumber(red->psubcomm,red->nsubcomm);CHKERRQ(ierr); 366d8a68f86SHong Zhang ierr = PetscSubcommSetType(red->psubcomm,PETSC_SUBCOMM_INTERLACED);CHKERRQ(ierr); 3675f06b7aaSBarry Smith ierr = PetscLogObjectMemory(pc,sizeof(PetscSubcomm));CHKERRQ(ierr); 3685f06b7aaSBarry Smith 3695f06b7aaSBarry Smith /* create a new PC that processors in each subcomm have copy of */ 3705f06b7aaSBarry Smith subcomm = red->psubcomm->comm; 3715f06b7aaSBarry Smith ierr = KSPCreate(subcomm,&red->ksp);CHKERRQ(ierr); 3725f06b7aaSBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)red->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 3735f06b7aaSBarry Smith ierr = PetscLogObjectParent(pc,red->ksp);CHKERRQ(ierr); 3745f06b7aaSBarry Smith ierr = KSPSetType(red->ksp,KSPPREONLY);CHKERRQ(ierr); 3755f06b7aaSBarry Smith ierr = KSPGetPC(red->ksp,&red->pc);CHKERRQ(ierr); 3765f06b7aaSBarry Smith ierr = PCSetType(red->pc,PCLU);CHKERRQ(ierr); 3775f06b7aaSBarry Smith 3785f06b7aaSBarry Smith ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 3795f06b7aaSBarry Smith ierr = KSPSetOptionsPrefix(red->ksp,prefix);CHKERRQ(ierr); 3805f06b7aaSBarry Smith ierr = KSPAppendOptionsPrefix(red->ksp,"redundant_");CHKERRQ(ierr); 3815f06b7aaSBarry Smith } 38283ab6a24SBarry Smith *innerksp = red->ksp; 3834b9ad928SBarry Smith PetscFunctionReturn(0); 3844b9ad928SBarry Smith } 3854b9ad928SBarry Smith EXTERN_C_END 3864b9ad928SBarry Smith 3874b9ad928SBarry Smith #undef __FUNCT__ 38883ab6a24SBarry Smith #define __FUNCT__ "PCRedundantGetKSP" 3894b9ad928SBarry Smith /*@ 39083ab6a24SBarry Smith PCRedundantGetKSP - Gets the less parallel KSP created by the redundant PC. 3914b9ad928SBarry Smith 3924b9ad928SBarry Smith Not Collective 3934b9ad928SBarry Smith 3944b9ad928SBarry Smith Input Parameter: 3954b9ad928SBarry Smith . pc - the preconditioner context 3964b9ad928SBarry Smith 3974b9ad928SBarry Smith Output Parameter: 39883ab6a24SBarry Smith . innerksp - the KSP on the smaller set of processes 3994b9ad928SBarry Smith 4004b9ad928SBarry Smith Level: advanced 4014b9ad928SBarry Smith 4024b9ad928SBarry Smith .keywords: PC, redundant solve 4034b9ad928SBarry Smith @*/ 40483ab6a24SBarry Smith PetscErrorCode PCRedundantGetKSP(PC pc,KSP *innerksp) 4054b9ad928SBarry Smith { 4064ac538c5SBarry Smith PetscErrorCode ierr; 4074b9ad928SBarry Smith 4084b9ad928SBarry Smith PetscFunctionBegin; 4090700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 41083ab6a24SBarry Smith PetscValidPointer(innerksp,2); 41183ab6a24SBarry Smith ierr = PetscTryMethod(pc,"PCRedundantGetKSP_C",(PC,KSP*),(pc,innerksp));CHKERRQ(ierr); 4124b9ad928SBarry Smith PetscFunctionReturn(0); 4134b9ad928SBarry Smith } 4144b9ad928SBarry Smith 4154b9ad928SBarry Smith EXTERN_C_BEGIN 4164b9ad928SBarry Smith #undef __FUNCT__ 4174b9ad928SBarry Smith #define __FUNCT__ "PCRedundantGetOperators_Redundant" 4187087cfbeSBarry Smith PetscErrorCode PCRedundantGetOperators_Redundant(PC pc,Mat *mat,Mat *pmat) 4194b9ad928SBarry Smith { 4204b9ad928SBarry Smith PC_Redundant *red = (PC_Redundant*)pc->data; 4214b9ad928SBarry Smith 4224b9ad928SBarry Smith PetscFunctionBegin; 423b3804887SHong Zhang if (mat) *mat = red->pmats; 424b3804887SHong Zhang if (pmat) *pmat = red->pmats; 4254b9ad928SBarry Smith PetscFunctionReturn(0); 4264b9ad928SBarry Smith } 4274b9ad928SBarry Smith EXTERN_C_END 4284b9ad928SBarry Smith 4294b9ad928SBarry Smith #undef __FUNCT__ 4304b9ad928SBarry Smith #define __FUNCT__ "PCRedundantGetOperators" 4314b9ad928SBarry Smith /*@ 4324b9ad928SBarry Smith PCRedundantGetOperators - gets the sequential matrix and preconditioner matrix 4334b9ad928SBarry Smith 4344b9ad928SBarry Smith Not Collective 4354b9ad928SBarry Smith 4364b9ad928SBarry Smith Input Parameter: 4374b9ad928SBarry Smith . pc - the preconditioner context 4384b9ad928SBarry Smith 4394b9ad928SBarry Smith Output Parameters: 4404b9ad928SBarry Smith + mat - the matrix 4414b9ad928SBarry Smith - pmat - the (possibly different) preconditioner matrix 4424b9ad928SBarry Smith 4434b9ad928SBarry Smith Level: advanced 4444b9ad928SBarry Smith 4454b9ad928SBarry Smith .keywords: PC, redundant solve 4464b9ad928SBarry Smith @*/ 4477087cfbeSBarry Smith PetscErrorCode PCRedundantGetOperators(PC pc,Mat *mat,Mat *pmat) 4484b9ad928SBarry Smith { 4494ac538c5SBarry Smith PetscErrorCode ierr; 4504b9ad928SBarry Smith 4514b9ad928SBarry Smith PetscFunctionBegin; 4520700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 4534482741eSBarry Smith if (mat) PetscValidPointer(mat,2); 4544482741eSBarry Smith if (pmat) PetscValidPointer(pmat,3); 4554ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCRedundantGetOperators_C",(PC,Mat*,Mat*),(pc,mat,pmat));CHKERRQ(ierr); 4564b9ad928SBarry Smith PetscFunctionReturn(0); 4574b9ad928SBarry Smith } 4584b9ad928SBarry Smith 4594b9ad928SBarry Smith /* -------------------------------------------------------------------------------------*/ 46037a17b4dSBarry Smith /*MC 46183ab6a24SBarry Smith PCREDUNDANT - Runs a KSP solver with preconditioner for the entire problem on subgroups of processors 46237a17b4dSBarry Smith 46383ab6a24SBarry Smith Options for the redundant preconditioners can be set with -redundant_pc_xxx for the redundant KSP with -redundant_ksp_xxx 46437a17b4dSBarry Smith 46509391456SBarry Smith Options Database: 4669b21d695SBarry Smith . -pc_redundant_number <n> - number of redundant solves, for example if you are using 64 MPI processes and 4679b21d695SBarry Smith use an n of 4 there will be 4 parallel solves each on 16 = 64/4 processes. 46809391456SBarry Smith 46937a17b4dSBarry Smith Level: intermediate 47037a17b4dSBarry Smith 47183ab6a24SBarry Smith Notes: The default KSP is preonly and the default PC is LU. 47283ab6a24SBarry Smith 47383ab6a24SBarry Smith Developer Notes: Note that PCSetInitialGuessNonzero() is not used by this class but likely should be. 4749cfaa89bSBarry Smith 47537a17b4dSBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PCRedundantSetScatter(), 47683ab6a24SBarry Smith PCRedundantGetKSP(), PCRedundantGetOperators(), PCRedundantSetNumber() 47737a17b4dSBarry Smith M*/ 47837a17b4dSBarry Smith 4794b9ad928SBarry Smith EXTERN_C_BEGIN 4804b9ad928SBarry Smith #undef __FUNCT__ 4814b9ad928SBarry Smith #define __FUNCT__ "PCCreate_Redundant" 4827087cfbeSBarry Smith PetscErrorCode PCCreate_Redundant(PC pc) 4834b9ad928SBarry Smith { 484dfbe8321SBarry Smith PetscErrorCode ierr; 4854b9ad928SBarry Smith PC_Redundant *red; 48669db28dcSHong Zhang PetscMPIInt size; 4873f457be1SHong Zhang 4884b9ad928SBarry Smith PetscFunctionBegin; 48938f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_Redundant,&red);CHKERRQ(ierr); 4907adad957SLisandro Dalcin ierr = MPI_Comm_size(((PetscObject)pc)->comm,&size);CHKERRQ(ierr); 49169db28dcSHong Zhang red->nsubcomm = size; 4924b9ad928SBarry Smith red->useparallelmat = PETSC_TRUE; 4931fbd8f88SHong Zhang pc->data = (void*)red; 4944b9ad928SBarry Smith 4954b9ad928SBarry Smith pc->ops->apply = PCApply_Redundant; 4964b9ad928SBarry Smith pc->ops->applytranspose = 0; 4974b9ad928SBarry Smith pc->ops->setup = PCSetUp_Redundant; 4984b9ad928SBarry Smith pc->ops->destroy = PCDestroy_Redundant; 4991ea5a559SBarry Smith pc->ops->reset = PCReset_Redundant; 5004b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_Redundant; 5014b9ad928SBarry Smith pc->ops->view = PCView_Redundant; 5024b9ad928SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCRedundantSetScatter_C","PCRedundantSetScatter_Redundant", 5034b9ad928SBarry Smith PCRedundantSetScatter_Redundant);CHKERRQ(ierr); 50409a6bc64SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCRedundantSetNumber_C","PCRedundantSetNumber_Redundant", 50509a6bc64SHong Zhang PCRedundantSetNumber_Redundant);CHKERRQ(ierr); 50683ab6a24SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCRedundantGetKSP_C","PCRedundantGetKSP_Redundant", 50783ab6a24SBarry Smith PCRedundantGetKSP_Redundant);CHKERRQ(ierr); 5084b9ad928SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCRedundantGetOperators_C","PCRedundantGetOperators_Redundant", 5094b9ad928SBarry Smith PCRedundantGetOperators_Redundant);CHKERRQ(ierr); 5104b9ad928SBarry Smith PetscFunctionReturn(0); 5114b9ad928SBarry Smith } 5124b9ad928SBarry Smith EXTERN_C_END 513