1dba47a55SKris Buschelman 24b9ad928SBarry Smith /* 33f457be1SHong Zhang This file defines a "solve the problem redundantly on each subgroup of processor" preconditioner. 44b9ad928SBarry Smith */ 507475bc1SBarry Smith #include <petsc-private/pcimpl.h> 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 */ 11ce94432eSBarry Smith Vec xsub,ysub; /* vectors of a subcommunicator to hold parallel vectors of PetscObjectComm((PetscObject)pc) */ 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; 30251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 31251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((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 } 484b9ad928SBarry Smith PetscFunctionReturn(0); 494b9ad928SBarry Smith } 504b9ad928SBarry Smith 51*d3b23db5SHong Zhang extern PetscErrorCode MatGetRedundantMatrix_MPIAIJ_interlaced(Mat,PetscInt,PetscSubcomm,MatReuse,Mat*); /* rm later! */ 52*d3b23db5SHong Zhang 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; 63ce94432eSBarry Smith MPI_Comm 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; 70ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 7123ce1328SBarry Smith ierr = MatGetVecs(pc->pmat,&vec,0);CHKERRQ(ierr); 7223ce1328SBarry Smith ierr = VecGetSize(vec,&m);CHKERRQ(ierr); 731fbd8f88SHong Zhang 744b9ad928SBarry Smith if (!pc->setupcalled) { 755f06b7aaSBarry Smith if (!red->psubcomm) { 76d8a68f86SHong Zhang ierr = PetscSubcommCreate(comm,&red->psubcomm);CHKERRQ(ierr); 77d8a68f86SHong Zhang ierr = PetscSubcommSetNumber(red->psubcomm,red->nsubcomm);CHKERRQ(ierr); 78d8a68f86SHong Zhang ierr = PetscSubcommSetType(red->psubcomm,PETSC_SUBCOMM_INTERLACED);CHKERRQ(ierr); 791fbd8f88SHong Zhang ierr = PetscLogObjectMemory(pc,sizeof(PetscSubcomm));CHKERRQ(ierr); 801fbd8f88SHong Zhang 811fbd8f88SHong Zhang /* create a new PC that processors in each subcomm have copy of */ 820d7810c8SBarry Smith subcomm = red->psubcomm->comm; 832fa5cd67SKarl Rupp 845f06b7aaSBarry Smith ierr = KSPCreate(subcomm,&red->ksp);CHKERRQ(ierr); 855f06b7aaSBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)red->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 865f06b7aaSBarry Smith ierr = PetscLogObjectParent(pc,red->ksp);CHKERRQ(ierr); 875f06b7aaSBarry Smith ierr = KSPSetType(red->ksp,KSPPREONLY);CHKERRQ(ierr); 885f06b7aaSBarry Smith ierr = KSPGetPC(red->ksp,&red->pc);CHKERRQ(ierr); 89cf52b8b1SHong Zhang ierr = PCSetType(red->pc,PCLU);CHKERRQ(ierr); 90cf52b8b1SHong Zhang 911fbd8f88SHong Zhang ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 925f06b7aaSBarry Smith ierr = KSPSetOptionsPrefix(red->ksp,prefix);CHKERRQ(ierr); 935f06b7aaSBarry Smith ierr = KSPAppendOptionsPrefix(red->ksp,"redundant_");CHKERRQ(ierr); 942fa5cd67SKarl Rupp } else subcomm = red->psubcomm->comm; 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 */ 1052fa5cd67SKarl Rupp if (subrank+1 < subsize) rend_sub = range[red->psubcomm->n*(subrank+1)]; 1062fa5cd67SKarl Rupp else rend_sub = m; 1072fa5cd67SKarl Rupp 1083f457be1SHong Zhang mloc_sub = rend_sub - rstart_sub; 1091fbd8f88SHong Zhang ierr = VecCreateMPI(subcomm,mloc_sub,PETSC_DECIDE,&red->ysub);CHKERRQ(ierr); 1103f457be1SHong Zhang /* create xsub with empty local arrays, because xdup's arrays will be placed into it */ 1110298fd71SBarry Smith ierr = VecCreateMPIWithArray(subcomm,1,mloc_sub,PETSC_DECIDE,NULL,&red->xsub);CHKERRQ(ierr); 1123f457be1SHong Zhang 1133f457be1SHong Zhang /* create xdup and ydup. ydup has empty local arrays because ysub's arrays will be place into it. 114ce94432eSBarry Smith Note: we use communicator dupcomm, not PetscObjectComm((PetscObject)pc)! */ 1151fbd8f88SHong Zhang ierr = VecCreateMPI(red->psubcomm->dupparent,mloc_sub,PETSC_DECIDE,&red->xdup);CHKERRQ(ierr); 1160298fd71SBarry Smith ierr = VecCreateMPIWithArray(red->psubcomm->dupparent,1,mloc_sub,PETSC_DECIDE,NULL,&red->ydup);CHKERRQ(ierr); 1173f457be1SHong Zhang 1183f457be1SHong Zhang /* create vec scatters */ 1193f457be1SHong Zhang if (!red->scatterin) { 1203f457be1SHong Zhang IS is1,is2; 1213f457be1SHong Zhang PetscInt *idx1,*idx2,i,j,k; 12245fc02eaSBarry Smith 1231d79065fSBarry Smith ierr = PetscMalloc2(red->psubcomm->n*mlocal,PetscInt,&idx1,red->psubcomm->n*mlocal,PetscInt,&idx2);CHKERRQ(ierr); 1243f457be1SHong Zhang j = 0; 1251fbd8f88SHong Zhang for (k=0; k<red->psubcomm->n; k++) { 1263f457be1SHong Zhang for (i=mstart; i<mend; i++) { 1273f457be1SHong Zhang idx1[j] = i; 1283f457be1SHong Zhang idx2[j++] = i + m*k; 1293f457be1SHong Zhang } 1303f457be1SHong Zhang } 13170b3c8c7SBarry Smith ierr = ISCreateGeneral(comm,red->psubcomm->n*mlocal,idx1,PETSC_COPY_VALUES,&is1);CHKERRQ(ierr); 13270b3c8c7SBarry Smith ierr = ISCreateGeneral(comm,red->psubcomm->n*mlocal,idx2,PETSC_COPY_VALUES,&is2);CHKERRQ(ierr); 1333f457be1SHong Zhang ierr = VecScatterCreate(vec,is1,red->xdup,is2,&red->scatterin);CHKERRQ(ierr); 134fcfd50ebSBarry Smith ierr = ISDestroy(&is1);CHKERRQ(ierr); 135fcfd50ebSBarry Smith ierr = ISDestroy(&is2);CHKERRQ(ierr); 1363f457be1SHong Zhang 1371fbd8f88SHong Zhang ierr = ISCreateStride(comm,mlocal,mstart+ red->psubcomm->color*m,1,&is1);CHKERRQ(ierr); 1383f457be1SHong Zhang ierr = ISCreateStride(comm,mlocal,mstart,1,&is2);CHKERRQ(ierr); 1393f457be1SHong Zhang ierr = VecScatterCreate(red->xdup,is1,vec,is2,&red->scatterout);CHKERRQ(ierr); 140fcfd50ebSBarry Smith ierr = ISDestroy(&is1);CHKERRQ(ierr); 141fcfd50ebSBarry Smith ierr = ISDestroy(&is2);CHKERRQ(ierr); 1421d79065fSBarry Smith ierr = PetscFree2(idx1,idx2);CHKERRQ(ierr); 1434b9ad928SBarry Smith } 1444b9ad928SBarry Smith } 1456bf464f9SBarry Smith ierr = VecDestroy(&vec);CHKERRQ(ierr); 1464b9ad928SBarry Smith 1474b9ad928SBarry Smith /* if pmatrix set by user is sequential then we do not need to gather the parallel matrix */ 1483f457be1SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1492fa5cd67SKarl Rupp if (size == 1) red->useparallelmat = PETSC_FALSE; 1504b9ad928SBarry Smith 1514b9ad928SBarry Smith if (red->useparallelmat) { 1524b9ad928SBarry Smith if (pc->setupcalled == 1 && pc->flag == DIFFERENT_NONZERO_PATTERN) { 1534b9ad928SBarry Smith /* destroy old matrices */ 1546bf464f9SBarry Smith ierr = MatDestroy(&red->pmats);CHKERRQ(ierr); 1554b9ad928SBarry Smith } else if (pc->setupcalled == 1) { 1564b9ad928SBarry Smith reuse = MAT_REUSE_MATRIX; 1574b9ad928SBarry Smith str = SAME_NONZERO_PATTERN; 1584b9ad928SBarry Smith } 1594b9ad928SBarry Smith 1603f457be1SHong Zhang /* grab the parallel matrix and put it into processors of a subcomminicator */ 161f664ae05SHong Zhang /*--------------------------------------------------------------------------*/ 162f664ae05SHong Zhang ierr = VecGetLocalSize(red->ysub,&mlocal_sub);CHKERRQ(ierr); 163*d3b23db5SHong Zhang ierr = MatGetRedundantMatrix_MPIAIJ_interlaced(pc->pmat,red->psubcomm->n,red->psubcomm,reuse,&red->pmats);CHKERRQ(ierr); 1643f457be1SHong Zhang /* tell PC of the subcommunicator its operators */ 16590f1c854SHong Zhang ierr = KSPSetOperators(red->ksp,red->pmats,red->pmats,str);CHKERRQ(ierr); 1664b9ad928SBarry Smith } else { 16790f1c854SHong Zhang ierr = KSPSetOperators(red->ksp,pc->mat,pc->pmat,pc->flag);CHKERRQ(ierr); 1684b9ad928SBarry Smith } 1690c24e6a1SHong Zhang if (pc->setfromoptionscalled) { 1703e065800SHong Zhang ierr = KSPSetFromOptions(red->ksp);CHKERRQ(ierr); 1710c24e6a1SHong Zhang } 1723e065800SHong Zhang ierr = KSPSetUp(red->ksp);CHKERRQ(ierr); 1734b9ad928SBarry Smith PetscFunctionReturn(0); 1744b9ad928SBarry Smith } 1754b9ad928SBarry Smith 1764b9ad928SBarry Smith #undef __FUNCT__ 1774b9ad928SBarry Smith #define __FUNCT__ "PCApply_Redundant" 1786849ba73SBarry Smith static PetscErrorCode PCApply_Redundant(PC pc,Vec x,Vec y) 1794b9ad928SBarry Smith { 1804b9ad928SBarry Smith PC_Redundant *red = (PC_Redundant*)pc->data; 181dfbe8321SBarry Smith PetscErrorCode ierr; 1823f457be1SHong Zhang PetscScalar *array; 1834b9ad928SBarry Smith 1844b9ad928SBarry Smith PetscFunctionBegin; 1853f457be1SHong Zhang /* scatter x to xdup */ 186ca9f406cSSatish Balay ierr = VecScatterBegin(red->scatterin,x,red->xdup,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 187ca9f406cSSatish Balay ierr = VecScatterEnd(red->scatterin,x,red->xdup,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1883f457be1SHong Zhang 1893f457be1SHong Zhang /* place xdup's local array into xsub */ 1903f457be1SHong Zhang ierr = VecGetArray(red->xdup,&array);CHKERRQ(ierr); 1913f457be1SHong Zhang ierr = VecPlaceArray(red->xsub,(const PetscScalar*)array);CHKERRQ(ierr); 1924b9ad928SBarry Smith 1934b9ad928SBarry Smith /* apply preconditioner on each processor */ 19483ab6a24SBarry Smith ierr = KSPSolve(red->ksp,red->xsub,red->ysub);CHKERRQ(ierr); 1953f457be1SHong Zhang ierr = VecResetArray(red->xsub);CHKERRQ(ierr); 1963f457be1SHong Zhang ierr = VecRestoreArray(red->xdup,&array);CHKERRQ(ierr); 1974b9ad928SBarry Smith 1983f457be1SHong Zhang /* place ysub's local array into ydup */ 1993f457be1SHong Zhang ierr = VecGetArray(red->ysub,&array);CHKERRQ(ierr); 2003f457be1SHong Zhang ierr = VecPlaceArray(red->ydup,(const PetscScalar*)array);CHKERRQ(ierr); 2013f457be1SHong Zhang 2023f457be1SHong Zhang /* scatter ydup to y */ 203ca9f406cSSatish Balay ierr = VecScatterBegin(red->scatterout,red->ydup,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 204ca9f406cSSatish Balay ierr = VecScatterEnd(red->scatterout,red->ydup,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 2053f457be1SHong Zhang ierr = VecResetArray(red->ydup);CHKERRQ(ierr); 2063f457be1SHong Zhang ierr = VecRestoreArray(red->ysub,&array);CHKERRQ(ierr); 2074b9ad928SBarry Smith PetscFunctionReturn(0); 2084b9ad928SBarry Smith } 2094b9ad928SBarry Smith 2104b9ad928SBarry Smith #undef __FUNCT__ 2111ea5a559SBarry Smith #define __FUNCT__ "PCReset_Redundant" 2121ea5a559SBarry Smith static PetscErrorCode PCReset_Redundant(PC pc) 2134b9ad928SBarry Smith { 2144b9ad928SBarry Smith PC_Redundant *red = (PC_Redundant*)pc->data; 215dfbe8321SBarry Smith PetscErrorCode ierr; 2164b9ad928SBarry Smith 2174b9ad928SBarry Smith PetscFunctionBegin; 2186bf464f9SBarry Smith ierr = VecScatterDestroy(&red->scatterin);CHKERRQ(ierr); 2196bf464f9SBarry Smith ierr = VecScatterDestroy(&red->scatterout);CHKERRQ(ierr); 2206bf464f9SBarry Smith ierr = VecDestroy(&red->ysub);CHKERRQ(ierr); 2216bf464f9SBarry Smith ierr = VecDestroy(&red->xsub);CHKERRQ(ierr); 2226bf464f9SBarry Smith ierr = VecDestroy(&red->xdup);CHKERRQ(ierr); 2236bf464f9SBarry Smith ierr = VecDestroy(&red->ydup);CHKERRQ(ierr); 2246bf464f9SBarry Smith ierr = MatDestroy(&red->pmats);CHKERRQ(ierr); 225546c30cbSLisandro Dalcin if (red->ksp) {ierr = KSPReset(red->ksp);CHKERRQ(ierr);} 2261ea5a559SBarry Smith PetscFunctionReturn(0); 2271ea5a559SBarry Smith } 2281ea5a559SBarry Smith 2291ea5a559SBarry Smith #undef __FUNCT__ 2301ea5a559SBarry Smith #define __FUNCT__ "PCDestroy_Redundant" 2311ea5a559SBarry Smith static PetscErrorCode PCDestroy_Redundant(PC pc) 2321ea5a559SBarry Smith { 2331ea5a559SBarry Smith PC_Redundant *red = (PC_Redundant*)pc->data; 2341ea5a559SBarry Smith PetscErrorCode ierr; 2351ea5a559SBarry Smith 2361ea5a559SBarry Smith PetscFunctionBegin; 2371ea5a559SBarry Smith ierr = PCReset_Redundant(pc);CHKERRQ(ierr); 2386bf464f9SBarry Smith ierr = KSPDestroy(&red->ksp);CHKERRQ(ierr); 2396bf464f9SBarry Smith ierr = PetscSubcommDestroy(&red->psubcomm);CHKERRQ(ierr); 240c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 2414b9ad928SBarry Smith PetscFunctionReturn(0); 2424b9ad928SBarry Smith } 2434b9ad928SBarry Smith 2444b9ad928SBarry Smith #undef __FUNCT__ 2454b9ad928SBarry Smith #define __FUNCT__ "PCSetFromOptions_Redundant" 2466849ba73SBarry Smith static PetscErrorCode PCSetFromOptions_Redundant(PC pc) 2474b9ad928SBarry Smith { 248a98ce0f4SHong Zhang PetscErrorCode ierr; 249a98ce0f4SHong Zhang PC_Redundant *red = (PC_Redundant*)pc->data; 250a98ce0f4SHong Zhang 2514b9ad928SBarry Smith PetscFunctionBegin; 252a98ce0f4SHong Zhang ierr = PetscOptionsHead("Redundant options");CHKERRQ(ierr); 25309a6bc64SHong Zhang ierr = PetscOptionsInt("-pc_redundant_number","Number of redundant pc","PCRedundantSetNumber",red->nsubcomm,&red->nsubcomm,0);CHKERRQ(ierr); 254a98ce0f4SHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 2554b9ad928SBarry Smith PetscFunctionReturn(0); 2564b9ad928SBarry Smith } 2574b9ad928SBarry Smith 2584b9ad928SBarry Smith #undef __FUNCT__ 25909a6bc64SHong Zhang #define __FUNCT__ "PCRedundantSetNumber_Redundant" 260f7a08781SBarry Smith static PetscErrorCode PCRedundantSetNumber_Redundant(PC pc,PetscInt nreds) 26109a6bc64SHong Zhang { 26209a6bc64SHong Zhang PC_Redundant *red = (PC_Redundant*)pc->data; 26309a6bc64SHong Zhang 26409a6bc64SHong Zhang PetscFunctionBegin; 26509a6bc64SHong Zhang red->nsubcomm = nreds; 26609a6bc64SHong Zhang PetscFunctionReturn(0); 26709a6bc64SHong Zhang } 26809a6bc64SHong Zhang 26909a6bc64SHong Zhang #undef __FUNCT__ 27009a6bc64SHong Zhang #define __FUNCT__ "PCRedundantSetNumber" 27109a6bc64SHong Zhang /*@ 27209a6bc64SHong Zhang PCRedundantSetNumber - Sets the number of redundant preconditioner contexts. 27309a6bc64SHong Zhang 2743f9fe445SBarry Smith Logically Collective on PC 27509a6bc64SHong Zhang 27609a6bc64SHong Zhang Input Parameters: 27709a6bc64SHong Zhang + pc - the preconditioner context 2789b21d695SBarry Smith - nredundant - number of redundant preconditioner contexts; for example if you are using 64 MPI processes and 2799b21d695SBarry Smith use an nredundant of 4 there will be 4 parallel solves each on 16 = 64/4 processes. 28009a6bc64SHong Zhang 28109a6bc64SHong Zhang Level: advanced 28209a6bc64SHong Zhang 28309a6bc64SHong Zhang .keywords: PC, redundant solve 28409a6bc64SHong Zhang @*/ 2857087cfbeSBarry Smith PetscErrorCode PCRedundantSetNumber(PC pc,PetscInt nredundant) 28609a6bc64SHong Zhang { 2874ac538c5SBarry Smith PetscErrorCode ierr; 28809a6bc64SHong Zhang 28909a6bc64SHong Zhang PetscFunctionBegin; 2900700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 291ce94432eSBarry Smith if (nredundant <= 0) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG, "num of redundant pc %D must be positive",nredundant); 2924ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCRedundantSetNumber_C",(PC,PetscInt),(pc,nredundant));CHKERRQ(ierr); 29309a6bc64SHong Zhang PetscFunctionReturn(0); 29409a6bc64SHong Zhang } 29509a6bc64SHong Zhang 29609a6bc64SHong Zhang #undef __FUNCT__ 2974b9ad928SBarry Smith #define __FUNCT__ "PCRedundantSetScatter_Redundant" 298f7a08781SBarry Smith static PetscErrorCode PCRedundantSetScatter_Redundant(PC pc,VecScatter in,VecScatter out) 2994b9ad928SBarry Smith { 3004b9ad928SBarry Smith PC_Redundant *red = (PC_Redundant*)pc->data; 301dfbe8321SBarry Smith PetscErrorCode ierr; 3024b9ad928SBarry Smith 3034b9ad928SBarry Smith PetscFunctionBegin; 3044b9ad928SBarry Smith ierr = PetscObjectReference((PetscObject)in);CHKERRQ(ierr); 3056bf464f9SBarry Smith ierr = VecScatterDestroy(&red->scatterin);CHKERRQ(ierr); 3062fa5cd67SKarl Rupp 307c3122656SLisandro Dalcin red->scatterin = in; 3082fa5cd67SKarl Rupp 3094b9ad928SBarry Smith ierr = PetscObjectReference((PetscObject)out);CHKERRQ(ierr); 3106bf464f9SBarry Smith ierr = VecScatterDestroy(&red->scatterout);CHKERRQ(ierr); 311c3122656SLisandro Dalcin red->scatterout = out; 3124b9ad928SBarry Smith PetscFunctionReturn(0); 3134b9ad928SBarry Smith } 3144b9ad928SBarry Smith 3154b9ad928SBarry Smith #undef __FUNCT__ 3164b9ad928SBarry Smith #define __FUNCT__ "PCRedundantSetScatter" 3174b9ad928SBarry Smith /*@ 3184b9ad928SBarry Smith PCRedundantSetScatter - Sets the scatter used to copy values into the 3194b9ad928SBarry Smith redundant local solve and the scatter to move them back into the global 3204b9ad928SBarry Smith vector. 3214b9ad928SBarry Smith 3223f9fe445SBarry Smith Logically Collective on PC 3234b9ad928SBarry Smith 3244b9ad928SBarry Smith Input Parameters: 3254b9ad928SBarry Smith + pc - the preconditioner context 3264b9ad928SBarry Smith . in - the scatter to move the values in 3274b9ad928SBarry Smith - out - the scatter to move them out 3284b9ad928SBarry Smith 3294b9ad928SBarry Smith Level: advanced 3304b9ad928SBarry Smith 3314b9ad928SBarry Smith .keywords: PC, redundant solve 3324b9ad928SBarry Smith @*/ 3337087cfbeSBarry Smith PetscErrorCode PCRedundantSetScatter(PC pc,VecScatter in,VecScatter out) 3344b9ad928SBarry Smith { 3354ac538c5SBarry Smith PetscErrorCode ierr; 3364b9ad928SBarry Smith 3374b9ad928SBarry Smith PetscFunctionBegin; 3380700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 3390700a824SBarry Smith PetscValidHeaderSpecific(in,VEC_SCATTER_CLASSID,2); 3400700a824SBarry Smith PetscValidHeaderSpecific(out,VEC_SCATTER_CLASSID,3); 3414ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCRedundantSetScatter_C",(PC,VecScatter,VecScatter),(pc,in,out));CHKERRQ(ierr); 3424b9ad928SBarry Smith PetscFunctionReturn(0); 3434b9ad928SBarry Smith } 3444b9ad928SBarry Smith 3454b9ad928SBarry Smith #undef __FUNCT__ 34683ab6a24SBarry Smith #define __FUNCT__ "PCRedundantGetKSP_Redundant" 347f7a08781SBarry Smith static PetscErrorCode PCRedundantGetKSP_Redundant(PC pc,KSP *innerksp) 3484b9ad928SBarry Smith { 3495f06b7aaSBarry Smith PetscErrorCode ierr; 3504b9ad928SBarry Smith PC_Redundant *red = (PC_Redundant*)pc->data; 3515f06b7aaSBarry Smith MPI_Comm comm,subcomm; 3525f06b7aaSBarry Smith const char *prefix; 3534b9ad928SBarry Smith 3544b9ad928SBarry Smith PetscFunctionBegin; 3555f06b7aaSBarry Smith if (!red->psubcomm) { 3565f06b7aaSBarry Smith ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 357d8a68f86SHong Zhang ierr = PetscSubcommCreate(comm,&red->psubcomm);CHKERRQ(ierr); 358d8a68f86SHong Zhang ierr = PetscSubcommSetNumber(red->psubcomm,red->nsubcomm);CHKERRQ(ierr); 359d8a68f86SHong Zhang ierr = PetscSubcommSetType(red->psubcomm,PETSC_SUBCOMM_INTERLACED);CHKERRQ(ierr); 3605f06b7aaSBarry Smith ierr = PetscLogObjectMemory(pc,sizeof(PetscSubcomm));CHKERRQ(ierr); 3615f06b7aaSBarry Smith 3625f06b7aaSBarry Smith /* create a new PC that processors in each subcomm have copy of */ 3635f06b7aaSBarry Smith subcomm = red->psubcomm->comm; 3642fa5cd67SKarl Rupp 3655f06b7aaSBarry Smith ierr = KSPCreate(subcomm,&red->ksp);CHKERRQ(ierr); 3665f06b7aaSBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)red->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 3675f06b7aaSBarry Smith ierr = PetscLogObjectParent(pc,red->ksp);CHKERRQ(ierr); 3685f06b7aaSBarry Smith ierr = KSPSetType(red->ksp,KSPPREONLY);CHKERRQ(ierr); 3695f06b7aaSBarry Smith ierr = KSPGetPC(red->ksp,&red->pc);CHKERRQ(ierr); 3705f06b7aaSBarry Smith ierr = PCSetType(red->pc,PCLU);CHKERRQ(ierr); 3715f06b7aaSBarry Smith 3725f06b7aaSBarry Smith ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 3735f06b7aaSBarry Smith ierr = KSPSetOptionsPrefix(red->ksp,prefix);CHKERRQ(ierr); 3745f06b7aaSBarry Smith ierr = KSPAppendOptionsPrefix(red->ksp,"redundant_");CHKERRQ(ierr); 3755f06b7aaSBarry Smith } 37683ab6a24SBarry Smith *innerksp = red->ksp; 3774b9ad928SBarry Smith PetscFunctionReturn(0); 3784b9ad928SBarry Smith } 3794b9ad928SBarry Smith 3804b9ad928SBarry Smith #undef __FUNCT__ 38183ab6a24SBarry Smith #define __FUNCT__ "PCRedundantGetKSP" 3824b9ad928SBarry Smith /*@ 38383ab6a24SBarry Smith PCRedundantGetKSP - Gets the less parallel KSP created by the redundant PC. 3844b9ad928SBarry Smith 3854b9ad928SBarry Smith Not Collective 3864b9ad928SBarry Smith 3874b9ad928SBarry Smith Input Parameter: 3884b9ad928SBarry Smith . pc - the preconditioner context 3894b9ad928SBarry Smith 3904b9ad928SBarry Smith Output Parameter: 39183ab6a24SBarry Smith . innerksp - the KSP on the smaller set of processes 3924b9ad928SBarry Smith 3934b9ad928SBarry Smith Level: advanced 3944b9ad928SBarry Smith 3954b9ad928SBarry Smith .keywords: PC, redundant solve 3964b9ad928SBarry Smith @*/ 39783ab6a24SBarry Smith PetscErrorCode PCRedundantGetKSP(PC pc,KSP *innerksp) 3984b9ad928SBarry Smith { 3994ac538c5SBarry Smith PetscErrorCode ierr; 4004b9ad928SBarry Smith 4014b9ad928SBarry Smith PetscFunctionBegin; 4020700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 40383ab6a24SBarry Smith PetscValidPointer(innerksp,2); 40483ab6a24SBarry Smith ierr = PetscTryMethod(pc,"PCRedundantGetKSP_C",(PC,KSP*),(pc,innerksp));CHKERRQ(ierr); 4054b9ad928SBarry Smith PetscFunctionReturn(0); 4064b9ad928SBarry Smith } 4074b9ad928SBarry Smith 4084b9ad928SBarry Smith #undef __FUNCT__ 4094b9ad928SBarry Smith #define __FUNCT__ "PCRedundantGetOperators_Redundant" 410f7a08781SBarry Smith static PetscErrorCode PCRedundantGetOperators_Redundant(PC pc,Mat *mat,Mat *pmat) 4114b9ad928SBarry Smith { 4124b9ad928SBarry Smith PC_Redundant *red = (PC_Redundant*)pc->data; 4134b9ad928SBarry Smith 4144b9ad928SBarry Smith PetscFunctionBegin; 415b3804887SHong Zhang if (mat) *mat = red->pmats; 416b3804887SHong Zhang if (pmat) *pmat = red->pmats; 4174b9ad928SBarry Smith PetscFunctionReturn(0); 4184b9ad928SBarry Smith } 4194b9ad928SBarry Smith 4204b9ad928SBarry Smith #undef __FUNCT__ 4214b9ad928SBarry Smith #define __FUNCT__ "PCRedundantGetOperators" 4224b9ad928SBarry Smith /*@ 4234b9ad928SBarry Smith PCRedundantGetOperators - gets the sequential matrix and preconditioner matrix 4244b9ad928SBarry Smith 4254b9ad928SBarry Smith Not Collective 4264b9ad928SBarry Smith 4274b9ad928SBarry Smith Input Parameter: 4284b9ad928SBarry Smith . pc - the preconditioner context 4294b9ad928SBarry Smith 4304b9ad928SBarry Smith Output Parameters: 4314b9ad928SBarry Smith + mat - the matrix 4324b9ad928SBarry Smith - pmat - the (possibly different) preconditioner matrix 4334b9ad928SBarry Smith 4344b9ad928SBarry Smith Level: advanced 4354b9ad928SBarry Smith 4364b9ad928SBarry Smith .keywords: PC, redundant solve 4374b9ad928SBarry Smith @*/ 4387087cfbeSBarry Smith PetscErrorCode PCRedundantGetOperators(PC pc,Mat *mat,Mat *pmat) 4394b9ad928SBarry Smith { 4404ac538c5SBarry Smith PetscErrorCode ierr; 4414b9ad928SBarry Smith 4424b9ad928SBarry Smith PetscFunctionBegin; 4430700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 4444482741eSBarry Smith if (mat) PetscValidPointer(mat,2); 4454482741eSBarry Smith if (pmat) PetscValidPointer(pmat,3); 4464ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCRedundantGetOperators_C",(PC,Mat*,Mat*),(pc,mat,pmat));CHKERRQ(ierr); 4474b9ad928SBarry Smith PetscFunctionReturn(0); 4484b9ad928SBarry Smith } 4494b9ad928SBarry Smith 4504b9ad928SBarry Smith /* -------------------------------------------------------------------------------------*/ 45137a17b4dSBarry Smith /*MC 45283ab6a24SBarry Smith PCREDUNDANT - Runs a KSP solver with preconditioner for the entire problem on subgroups of processors 45337a17b4dSBarry Smith 45483ab6a24SBarry Smith Options for the redundant preconditioners can be set with -redundant_pc_xxx for the redundant KSP with -redundant_ksp_xxx 45537a17b4dSBarry Smith 45609391456SBarry Smith Options Database: 4579b21d695SBarry Smith . -pc_redundant_number <n> - number of redundant solves, for example if you are using 64 MPI processes and 4589b21d695SBarry Smith use an n of 4 there will be 4 parallel solves each on 16 = 64/4 processes. 45909391456SBarry Smith 46037a17b4dSBarry Smith Level: intermediate 46137a17b4dSBarry Smith 46283ab6a24SBarry Smith Notes: The default KSP is preonly and the default PC is LU. 46383ab6a24SBarry Smith 46483ab6a24SBarry Smith Developer Notes: Note that PCSetInitialGuessNonzero() is not used by this class but likely should be. 4659cfaa89bSBarry Smith 46637a17b4dSBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PCRedundantSetScatter(), 46783ab6a24SBarry Smith PCRedundantGetKSP(), PCRedundantGetOperators(), PCRedundantSetNumber() 46837a17b4dSBarry Smith M*/ 46937a17b4dSBarry Smith 4704b9ad928SBarry Smith #undef __FUNCT__ 4714b9ad928SBarry Smith #define __FUNCT__ "PCCreate_Redundant" 4728cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_Redundant(PC pc) 4734b9ad928SBarry Smith { 474dfbe8321SBarry Smith PetscErrorCode ierr; 4754b9ad928SBarry Smith PC_Redundant *red; 47669db28dcSHong Zhang PetscMPIInt size; 4773f457be1SHong Zhang 4784b9ad928SBarry Smith PetscFunctionBegin; 47938f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_Redundant,&red);CHKERRQ(ierr); 480ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr); 4812fa5cd67SKarl Rupp 48269db28dcSHong Zhang red->nsubcomm = size; 4834b9ad928SBarry Smith red->useparallelmat = PETSC_TRUE; 4841fbd8f88SHong Zhang pc->data = (void*)red; 4854b9ad928SBarry Smith 4864b9ad928SBarry Smith pc->ops->apply = PCApply_Redundant; 4874b9ad928SBarry Smith pc->ops->applytranspose = 0; 4884b9ad928SBarry Smith pc->ops->setup = PCSetUp_Redundant; 4894b9ad928SBarry Smith pc->ops->destroy = PCDestroy_Redundant; 4901ea5a559SBarry Smith pc->ops->reset = PCReset_Redundant; 4914b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_Redundant; 4924b9ad928SBarry Smith pc->ops->view = PCView_Redundant; 4932fa5cd67SKarl Rupp 494bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCRedundantSetScatter_C",PCRedundantSetScatter_Redundant);CHKERRQ(ierr); 495bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCRedundantSetNumber_C",PCRedundantSetNumber_Redundant);CHKERRQ(ierr); 496bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCRedundantGetKSP_C",PCRedundantGetKSP_Redundant);CHKERRQ(ierr); 497bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCRedundantGetOperators_C",PCRedundantGetOperators_Redundant);CHKERRQ(ierr); 4984b9ad928SBarry Smith PetscFunctionReturn(0); 4994b9ad928SBarry Smith } 500b2573a8aSBarry Smith 501