1dba47a55SKris Buschelman #define PETSCKSP_DLL 2dba47a55SKris Buschelman 34b9ad928SBarry Smith /* 43f457be1SHong Zhang This file defines a "solve the problem redundantly on each subgroup of processor" preconditioner. 54b9ad928SBarry Smith */ 66356e834SBarry Smith #include "private/pcimpl.h" /*I "petscpc.h" I*/ 74b9ad928SBarry Smith #include "petscksp.h" 84b9ad928SBarry Smith 94b9ad928SBarry Smith typedef struct { 103e065800SHong Zhang KSP ksp; 114b9ad928SBarry Smith PC pc; /* actual preconditioner used on each processor */ 127adad957SLisandro Dalcin Vec xsub,ysub; /* vectors of a subcommunicator to hold parallel vectors of ((PetscObject)pc)->comm */ 133f457be1SHong Zhang Vec xdup,ydup; /* parallel vector that congregates xsub or ysub facilitating vector scattering */ 14b3804887SHong Zhang Mat pmats; /* matrix and optional preconditioner matrix belong to a subcommunicator */ 153f457be1SHong Zhang VecScatter scatterin,scatterout; /* scatter used to move all values to each processor group (subcommunicator) */ 164b9ad928SBarry Smith PetscTruth useparallelmat; 17c540e29cSHong Zhang PetscSubcomm psubcomm; 181fbd8f88SHong Zhang PetscInt nsubcomm; /* num of data structure PetscSubcomm */ 194b9ad928SBarry Smith } PC_Redundant; 204b9ad928SBarry Smith 214b9ad928SBarry Smith #undef __FUNCT__ 224b9ad928SBarry Smith #define __FUNCT__ "PCView_Redundant" 236849ba73SBarry Smith static PetscErrorCode PCView_Redundant(PC pc,PetscViewer viewer) 244b9ad928SBarry Smith { 254b9ad928SBarry Smith PC_Redundant *red = (PC_Redundant*)pc->data; 26dfbe8321SBarry Smith PetscErrorCode ierr; 2732077d6dSBarry Smith PetscTruth iascii,isstring; 2803ccd0b4SBarry Smith PetscViewer subviewer; 294b9ad928SBarry Smith 304b9ad928SBarry Smith PetscFunctionBegin; 312692d6eeSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 322692d6eeSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSCVIEWERSTRING,&isstring);CHKERRQ(ierr); 3332077d6dSBarry Smith if (iascii) { 3403ccd0b4SBarry Smith if (!red->psubcomm) { 3503ccd0b4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Redundant preconditioner: Not yet setup\n");CHKERRQ(ierr); 3603ccd0b4SBarry Smith } else { 373e065800SHong Zhang ierr = PetscViewerASCIIPrintf(viewer," Redundant preconditioner: First (color=0) of %D PCs follows\n",red->nsubcomm);CHKERRQ(ierr); 387adad957SLisandro Dalcin ierr = PetscViewerGetSubcomm(viewer,((PetscObject)red->pc)->comm,&subviewer);CHKERRQ(ierr); 39*f5dd71faSBarry Smith if (!red->psubcomm->color) { /* only view first redundant pc */ 404b9ad928SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 413e065800SHong Zhang ierr = KSPView(red->ksp,subviewer);CHKERRQ(ierr); 424b9ad928SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 434b9ad928SBarry Smith } 447adad957SLisandro Dalcin ierr = PetscViewerRestoreSubcomm(viewer,((PetscObject)red->pc)->comm,&subviewer);CHKERRQ(ierr); 454b9ad928SBarry Smith } 4603ccd0b4SBarry Smith } else if (isstring) { 4703ccd0b4SBarry Smith ierr = PetscViewerStringSPrintf(viewer," Redundant solver preconditioner");CHKERRQ(ierr); 484b9ad928SBarry Smith } else { 4965e19b50SBarry Smith SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for PC redundant",((PetscObject)viewer)->type_name); 504b9ad928SBarry Smith } 514b9ad928SBarry Smith PetscFunctionReturn(0); 524b9ad928SBarry Smith } 534b9ad928SBarry Smith 544b9ad928SBarry Smith #undef __FUNCT__ 554b9ad928SBarry Smith #define __FUNCT__ "PCSetUp_Redundant" 566849ba73SBarry Smith static PetscErrorCode PCSetUp_Redundant(PC pc) 574b9ad928SBarry Smith { 584b9ad928SBarry Smith PC_Redundant *red = (PC_Redundant*)pc->data; 59dfbe8321SBarry Smith PetscErrorCode ierr; 6045fc02eaSBarry Smith PetscInt mstart,mend,mlocal,m,mlocal_sub,rstart_sub,rend_sub,mloc_sub; 6113f74950SBarry Smith PetscMPIInt size; 624b9ad928SBarry Smith MatReuse reuse = MAT_INITIAL_MATRIX; 634b9ad928SBarry Smith MatStructure str = DIFFERENT_NONZERO_PATTERN; 647adad957SLisandro Dalcin MPI_Comm comm = ((PetscObject)pc)->comm,subcomm; 6523ce1328SBarry Smith Vec vec; 663f457be1SHong Zhang PetscMPIInt subsize,subrank; 671fbd8f88SHong Zhang const char *prefix; 68b862ddfaSBarry Smith const PetscInt *range; 693f457be1SHong Zhang 704b9ad928SBarry Smith PetscFunctionBegin; 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) { 76638faf0bSHong Zhang ierr = PetscSubcommCreate(comm,red->nsubcomm,PETSC_SUBCOMM_INTERLACED,&red->psubcomm);CHKERRQ(ierr); 771fbd8f88SHong Zhang ierr = PetscLogObjectMemory(pc,sizeof(PetscSubcomm));CHKERRQ(ierr); 781fbd8f88SHong Zhang 791fbd8f88SHong Zhang /* create a new PC that processors in each subcomm have copy of */ 800d7810c8SBarry Smith subcomm = red->psubcomm->comm; 815f06b7aaSBarry Smith ierr = KSPCreate(subcomm,&red->ksp);CHKERRQ(ierr); 825f06b7aaSBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)red->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 835f06b7aaSBarry Smith ierr = PetscLogObjectParent(pc,red->ksp);CHKERRQ(ierr); 845f06b7aaSBarry Smith ierr = KSPSetType(red->ksp,KSPPREONLY);CHKERRQ(ierr); 855f06b7aaSBarry Smith ierr = KSPGetPC(red->ksp,&red->pc);CHKERRQ(ierr); 86cf52b8b1SHong Zhang ierr = PCSetType(red->pc,PCLU);CHKERRQ(ierr); 87cf52b8b1SHong Zhang 881fbd8f88SHong Zhang ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 895f06b7aaSBarry Smith ierr = KSPSetOptionsPrefix(red->ksp,prefix);CHKERRQ(ierr); 905f06b7aaSBarry Smith ierr = KSPAppendOptionsPrefix(red->ksp,"redundant_");CHKERRQ(ierr); 915f06b7aaSBarry Smith } else { 925f06b7aaSBarry Smith subcomm = red->psubcomm->comm; 935f06b7aaSBarry Smith } 941fbd8f88SHong Zhang 953f457be1SHong Zhang /* create working vectors xsub/ysub and xdup/ydup */ 9623ce1328SBarry Smith ierr = VecGetLocalSize(vec,&mlocal);CHKERRQ(ierr); 973f457be1SHong Zhang ierr = VecGetOwnershipRange(vec,&mstart,&mend);CHKERRQ(ierr); 984b9ad928SBarry Smith 993f457be1SHong Zhang /* get local size of xsub/ysub */ 1001fbd8f88SHong Zhang ierr = MPI_Comm_size(subcomm,&subsize);CHKERRQ(ierr); 1011fbd8f88SHong Zhang ierr = MPI_Comm_rank(subcomm,&subrank);CHKERRQ(ierr); 102b862ddfaSBarry Smith ierr = MatGetOwnershipRanges(pc->pmat,&range);CHKERRQ(ierr); 103b862ddfaSBarry Smith rstart_sub = range[red->psubcomm->n*subrank]; /* rstart in xsub/ysub */ 1043f457be1SHong Zhang if (subrank+1 < subsize){ 105b862ddfaSBarry Smith rend_sub = range[red->psubcomm->n*(subrank+1)]; 1063f457be1SHong Zhang } else { 1073f457be1SHong Zhang rend_sub = m; 1083f457be1SHong Zhang } 1093f457be1SHong Zhang mloc_sub = rend_sub - rstart_sub; 1101fbd8f88SHong Zhang ierr = VecCreateMPI(subcomm,mloc_sub,PETSC_DECIDE,&red->ysub);CHKERRQ(ierr); 1113f457be1SHong Zhang /* create xsub with empty local arrays, because xdup's arrays will be placed into it */ 1121fbd8f88SHong Zhang ierr = VecCreateMPIWithArray(subcomm,mloc_sub,PETSC_DECIDE,PETSC_NULL,&red->xsub);CHKERRQ(ierr); 1133f457be1SHong Zhang 1143f457be1SHong Zhang /* create xdup and ydup. ydup has empty local arrays because ysub's arrays will be place into it. 1157adad957SLisandro Dalcin Note: we use communicator dupcomm, not ((PetscObject)pc)->comm! */ 1161fbd8f88SHong Zhang ierr = VecCreateMPI(red->psubcomm->dupparent,mloc_sub,PETSC_DECIDE,&red->xdup);CHKERRQ(ierr); 1171fbd8f88SHong Zhang ierr = VecCreateMPIWithArray(red->psubcomm->dupparent,mloc_sub,PETSC_DECIDE,PETSC_NULL,&red->ydup);CHKERRQ(ierr); 1183f457be1SHong Zhang 1193f457be1SHong Zhang /* create vec scatters */ 1203f457be1SHong Zhang if (!red->scatterin){ 1213f457be1SHong Zhang IS is1,is2; 1223f457be1SHong Zhang PetscInt *idx1,*idx2,i,j,k; 12345fc02eaSBarry Smith 1241d79065fSBarry Smith ierr = PetscMalloc2(red->psubcomm->n*mlocal,PetscInt,&idx1,red->psubcomm->n*mlocal,PetscInt,&idx2);CHKERRQ(ierr); 1253f457be1SHong Zhang j = 0; 1261fbd8f88SHong Zhang for (k=0; k<red->psubcomm->n; k++){ 1273f457be1SHong Zhang for (i=mstart; i<mend; i++){ 1283f457be1SHong Zhang idx1[j] = i; 1293f457be1SHong Zhang idx2[j++] = i + m*k; 1303f457be1SHong Zhang } 1313f457be1SHong Zhang } 1321fbd8f88SHong Zhang ierr = ISCreateGeneral(comm,red->psubcomm->n*mlocal,idx1,&is1);CHKERRQ(ierr); 1331fbd8f88SHong Zhang ierr = ISCreateGeneral(comm,red->psubcomm->n*mlocal,idx2,&is2);CHKERRQ(ierr); 1343f457be1SHong Zhang ierr = VecScatterCreate(vec,is1,red->xdup,is2,&red->scatterin);CHKERRQ(ierr); 1353f457be1SHong Zhang ierr = ISDestroy(is1);CHKERRQ(ierr); 1363f457be1SHong Zhang ierr = ISDestroy(is2);CHKERRQ(ierr); 1373f457be1SHong Zhang 1381fbd8f88SHong Zhang ierr = ISCreateStride(comm,mlocal,mstart+ red->psubcomm->color*m,1,&is1);CHKERRQ(ierr); 1393f457be1SHong Zhang ierr = ISCreateStride(comm,mlocal,mstart,1,&is2);CHKERRQ(ierr); 1403f457be1SHong Zhang ierr = VecScatterCreate(red->xdup,is1,vec,is2,&red->scatterout);CHKERRQ(ierr); 1413f457be1SHong Zhang ierr = ISDestroy(is1);CHKERRQ(ierr); 1423f457be1SHong Zhang ierr = ISDestroy(is2);CHKERRQ(ierr); 1431d79065fSBarry Smith ierr = PetscFree2(idx1,idx2);CHKERRQ(ierr); 1444b9ad928SBarry Smith } 1454b9ad928SBarry Smith } 14623ce1328SBarry Smith ierr = VecDestroy(vec);CHKERRQ(ierr); 1474b9ad928SBarry Smith 1484b9ad928SBarry Smith /* if pmatrix set by user is sequential then we do not need to gather the parallel matrix */ 1493f457be1SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1504b9ad928SBarry Smith if (size == 1) { 1514b9ad928SBarry Smith red->useparallelmat = PETSC_FALSE; 1524b9ad928SBarry Smith } 1534b9ad928SBarry Smith 1544b9ad928SBarry Smith if (red->useparallelmat) { 1554b9ad928SBarry Smith if (pc->setupcalled == 1 && pc->flag == DIFFERENT_NONZERO_PATTERN) { 1564b9ad928SBarry Smith /* destroy old matrices */ 1574b9ad928SBarry Smith if (red->pmats) { 158b3804887SHong Zhang ierr = MatDestroy(red->pmats);CHKERRQ(ierr); 1594b9ad928SBarry Smith } 1604b9ad928SBarry Smith } else if (pc->setupcalled == 1) { 1614b9ad928SBarry Smith reuse = MAT_REUSE_MATRIX; 1624b9ad928SBarry Smith str = SAME_NONZERO_PATTERN; 1634b9ad928SBarry Smith } 1644b9ad928SBarry Smith 1653f457be1SHong Zhang /* grab the parallel matrix and put it into processors of a subcomminicator */ 166f664ae05SHong Zhang /*--------------------------------------------------------------------------*/ 167f664ae05SHong Zhang ierr = VecGetLocalSize(red->ysub,&mlocal_sub);CHKERRQ(ierr); 16869db28dcSHong Zhang ierr = MatGetRedundantMatrix(pc->pmat,red->psubcomm->n,red->psubcomm->comm,mlocal_sub,reuse,&red->pmats);CHKERRQ(ierr); 1693f457be1SHong Zhang /* tell PC of the subcommunicator its operators */ 17090f1c854SHong Zhang ierr = KSPSetOperators(red->ksp,red->pmats,red->pmats,str);CHKERRQ(ierr); 1714b9ad928SBarry Smith } else { 17290f1c854SHong Zhang ierr = KSPSetOperators(red->ksp,pc->mat,pc->pmat,pc->flag);CHKERRQ(ierr); 1734b9ad928SBarry Smith } 1740c24e6a1SHong Zhang if (pc->setfromoptionscalled){ 1753e065800SHong Zhang ierr = KSPSetFromOptions(red->ksp);CHKERRQ(ierr); 1760c24e6a1SHong Zhang } 1773e065800SHong Zhang ierr = KSPSetUp(red->ksp);CHKERRQ(ierr); 1784b9ad928SBarry Smith PetscFunctionReturn(0); 1794b9ad928SBarry Smith } 1804b9ad928SBarry Smith 1814b9ad928SBarry Smith #undef __FUNCT__ 1824b9ad928SBarry Smith #define __FUNCT__ "PCApply_Redundant" 1836849ba73SBarry Smith static PetscErrorCode PCApply_Redundant(PC pc,Vec x,Vec y) 1844b9ad928SBarry Smith { 1854b9ad928SBarry Smith PC_Redundant *red = (PC_Redundant*)pc->data; 186dfbe8321SBarry Smith PetscErrorCode ierr; 1873f457be1SHong Zhang PetscScalar *array; 1884b9ad928SBarry Smith 1894b9ad928SBarry Smith PetscFunctionBegin; 1903f457be1SHong Zhang /* scatter x to xdup */ 191ca9f406cSSatish Balay ierr = VecScatterBegin(red->scatterin,x,red->xdup,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 192ca9f406cSSatish Balay ierr = VecScatterEnd(red->scatterin,x,red->xdup,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1933f457be1SHong Zhang 1943f457be1SHong Zhang /* place xdup's local array into xsub */ 1953f457be1SHong Zhang ierr = VecGetArray(red->xdup,&array);CHKERRQ(ierr); 1963f457be1SHong Zhang ierr = VecPlaceArray(red->xsub,(const PetscScalar*)array);CHKERRQ(ierr); 1974b9ad928SBarry Smith 1984b9ad928SBarry Smith /* apply preconditioner on each processor */ 1993f457be1SHong Zhang ierr = PCApply(red->pc,red->xsub,red->ysub);CHKERRQ(ierr); 2003f457be1SHong Zhang ierr = VecResetArray(red->xsub);CHKERRQ(ierr); 2013f457be1SHong Zhang ierr = VecRestoreArray(red->xdup,&array);CHKERRQ(ierr); 2024b9ad928SBarry Smith 2033f457be1SHong Zhang /* place ysub's local array into ydup */ 2043f457be1SHong Zhang ierr = VecGetArray(red->ysub,&array);CHKERRQ(ierr); 2053f457be1SHong Zhang ierr = VecPlaceArray(red->ydup,(const PetscScalar*)array);CHKERRQ(ierr); 2063f457be1SHong Zhang 2073f457be1SHong Zhang /* scatter ydup to y */ 208ca9f406cSSatish Balay ierr = VecScatterBegin(red->scatterout,red->ydup,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 209ca9f406cSSatish Balay ierr = VecScatterEnd(red->scatterout,red->ydup,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 2103f457be1SHong Zhang ierr = VecResetArray(red->ydup);CHKERRQ(ierr); 2113f457be1SHong Zhang ierr = VecRestoreArray(red->ysub,&array);CHKERRQ(ierr); 2124b9ad928SBarry Smith PetscFunctionReturn(0); 2134b9ad928SBarry Smith } 2144b9ad928SBarry Smith 2154b9ad928SBarry Smith #undef __FUNCT__ 2164b9ad928SBarry Smith #define __FUNCT__ "PCDestroy_Redundant" 2176849ba73SBarry Smith static PetscErrorCode PCDestroy_Redundant(PC pc) 2184b9ad928SBarry Smith { 2194b9ad928SBarry Smith PC_Redundant *red = (PC_Redundant*)pc->data; 220dfbe8321SBarry Smith PetscErrorCode ierr; 2214b9ad928SBarry Smith 2224b9ad928SBarry Smith PetscFunctionBegin; 2234b9ad928SBarry Smith if (red->scatterin) {ierr = VecScatterDestroy(red->scatterin);CHKERRQ(ierr);} 2244b9ad928SBarry Smith if (red->scatterout) {ierr = VecScatterDestroy(red->scatterout);CHKERRQ(ierr);} 2253f457be1SHong Zhang if (red->ysub) {ierr = VecDestroy(red->ysub);CHKERRQ(ierr);} 2263f457be1SHong Zhang if (red->xsub) {ierr = VecDestroy(red->xsub);CHKERRQ(ierr);} 2273f457be1SHong Zhang if (red->xdup) {ierr = VecDestroy(red->xdup);CHKERRQ(ierr);} 2283f457be1SHong Zhang if (red->ydup) {ierr = VecDestroy(red->ydup);CHKERRQ(ierr);} 229b3804887SHong Zhang if (red->pmats) { 230b3804887SHong Zhang ierr = MatDestroy(red->pmats);CHKERRQ(ierr); 2313f457be1SHong Zhang } 2323e065800SHong Zhang if (red->ksp) {ierr = KSPDestroy(red->ksp);CHKERRQ(ierr);} 23330ca954eSBarry Smith if (red->psubcomm) {ierr = PetscSubcommDestroy(red->psubcomm);CHKERRQ(ierr);} 2344b9ad928SBarry Smith ierr = PetscFree(red);CHKERRQ(ierr); 2354b9ad928SBarry Smith PetscFunctionReturn(0); 2364b9ad928SBarry Smith } 2374b9ad928SBarry Smith 2384b9ad928SBarry Smith #undef __FUNCT__ 2394b9ad928SBarry Smith #define __FUNCT__ "PCSetFromOptions_Redundant" 2406849ba73SBarry Smith static PetscErrorCode PCSetFromOptions_Redundant(PC pc) 2414b9ad928SBarry Smith { 242a98ce0f4SHong Zhang PetscErrorCode ierr; 243a98ce0f4SHong Zhang PC_Redundant *red = (PC_Redundant*)pc->data; 244a98ce0f4SHong Zhang 2454b9ad928SBarry Smith PetscFunctionBegin; 246a98ce0f4SHong Zhang ierr = PetscOptionsHead("Redundant options");CHKERRQ(ierr); 24709a6bc64SHong Zhang ierr = PetscOptionsInt("-pc_redundant_number","Number of redundant pc","PCRedundantSetNumber",red->nsubcomm,&red->nsubcomm,0);CHKERRQ(ierr); 248a98ce0f4SHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 2494b9ad928SBarry Smith PetscFunctionReturn(0); 2504b9ad928SBarry Smith } 2514b9ad928SBarry Smith 2524b9ad928SBarry Smith EXTERN_C_BEGIN 2534b9ad928SBarry Smith #undef __FUNCT__ 25409a6bc64SHong Zhang #define __FUNCT__ "PCRedundantSetNumber_Redundant" 25509a6bc64SHong Zhang PetscErrorCode PETSCKSP_DLLEXPORT PCRedundantSetNumber_Redundant(PC pc,PetscInt nreds) 25609a6bc64SHong Zhang { 25709a6bc64SHong Zhang PC_Redundant *red = (PC_Redundant*)pc->data; 25809a6bc64SHong Zhang 25909a6bc64SHong Zhang PetscFunctionBegin; 26009a6bc64SHong Zhang red->nsubcomm = nreds; 26109a6bc64SHong Zhang PetscFunctionReturn(0); 26209a6bc64SHong Zhang } 26309a6bc64SHong Zhang EXTERN_C_END 26409a6bc64SHong Zhang 26509a6bc64SHong Zhang #undef __FUNCT__ 26609a6bc64SHong Zhang #define __FUNCT__ "PCRedundantSetNumber" 26709a6bc64SHong Zhang /*@ 26809a6bc64SHong Zhang PCRedundantSetNumber - Sets the number of redundant preconditioner contexts. 26909a6bc64SHong Zhang 2703f9fe445SBarry Smith Logically Collective on PC 27109a6bc64SHong Zhang 27209a6bc64SHong Zhang Input Parameters: 27309a6bc64SHong Zhang + pc - the preconditioner context 2749b21d695SBarry Smith - nredundant - number of redundant preconditioner contexts; for example if you are using 64 MPI processes and 2759b21d695SBarry Smith use an nredundant of 4 there will be 4 parallel solves each on 16 = 64/4 processes. 27609a6bc64SHong Zhang 27709a6bc64SHong Zhang Level: advanced 27809a6bc64SHong Zhang 27909a6bc64SHong Zhang .keywords: PC, redundant solve 28009a6bc64SHong Zhang @*/ 28109a6bc64SHong Zhang PetscErrorCode PETSCKSP_DLLEXPORT PCRedundantSetNumber(PC pc,PetscInt nredundant) 28209a6bc64SHong Zhang { 28309a6bc64SHong Zhang PetscErrorCode ierr,(*f)(PC,PetscInt); 28409a6bc64SHong Zhang 28509a6bc64SHong Zhang PetscFunctionBegin; 2860700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 28765e19b50SBarry Smith if (nredundant <= 0) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG, "num of redundant pc %D must be positive",nredundant); 28809a6bc64SHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCRedundantSetNumber_C",(void (**)(void))&f);CHKERRQ(ierr); 28909a6bc64SHong Zhang if (f) { 29009a6bc64SHong Zhang ierr = (*f)(pc,nredundant);CHKERRQ(ierr); 29109a6bc64SHong Zhang } 29209a6bc64SHong Zhang PetscFunctionReturn(0); 29309a6bc64SHong Zhang } 29409a6bc64SHong Zhang 29509a6bc64SHong Zhang EXTERN_C_BEGIN 29609a6bc64SHong Zhang #undef __FUNCT__ 2974b9ad928SBarry Smith #define __FUNCT__ "PCRedundantSetScatter_Redundant" 298dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT 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); 305c3122656SLisandro Dalcin if (red->scatterin) { ierr = VecScatterDestroy(red->scatterin);CHKERRQ(ierr); } 306c3122656SLisandro Dalcin red->scatterin = in; 3074b9ad928SBarry Smith ierr = PetscObjectReference((PetscObject)out);CHKERRQ(ierr); 308c3122656SLisandro Dalcin if (red->scatterout) { ierr = VecScatterDestroy(red->scatterout);CHKERRQ(ierr); } 309c3122656SLisandro Dalcin red->scatterout = out; 3104b9ad928SBarry Smith PetscFunctionReturn(0); 3114b9ad928SBarry Smith } 3124b9ad928SBarry Smith EXTERN_C_END 3134b9ad928SBarry Smith 3144b9ad928SBarry Smith #undef __FUNCT__ 3154b9ad928SBarry Smith #define __FUNCT__ "PCRedundantSetScatter" 3164b9ad928SBarry Smith /*@ 3174b9ad928SBarry Smith PCRedundantSetScatter - Sets the scatter used to copy values into the 3184b9ad928SBarry Smith redundant local solve and the scatter to move them back into the global 3194b9ad928SBarry Smith vector. 3204b9ad928SBarry Smith 3213f9fe445SBarry Smith Logically Collective on PC 3224b9ad928SBarry Smith 3234b9ad928SBarry Smith Input Parameters: 3244b9ad928SBarry Smith + pc - the preconditioner context 3254b9ad928SBarry Smith . in - the scatter to move the values in 3264b9ad928SBarry Smith - out - the scatter to move them out 3274b9ad928SBarry Smith 3284b9ad928SBarry Smith Level: advanced 3294b9ad928SBarry Smith 3304b9ad928SBarry Smith .keywords: PC, redundant solve 3314b9ad928SBarry Smith @*/ 332dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCRedundantSetScatter(PC pc,VecScatter in,VecScatter out) 3334b9ad928SBarry Smith { 334dfbe8321SBarry Smith PetscErrorCode ierr,(*f)(PC,VecScatter,VecScatter); 3354b9ad928SBarry Smith 3364b9ad928SBarry Smith PetscFunctionBegin; 3370700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 3380700a824SBarry Smith PetscValidHeaderSpecific(in,VEC_SCATTER_CLASSID,2); 3390700a824SBarry Smith PetscValidHeaderSpecific(out,VEC_SCATTER_CLASSID,3); 3404b9ad928SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCRedundantSetScatter_C",(void (**)(void))&f);CHKERRQ(ierr); 3414b9ad928SBarry Smith if (f) { 3424b9ad928SBarry Smith ierr = (*f)(pc,in,out);CHKERRQ(ierr); 3434b9ad928SBarry Smith } 3444b9ad928SBarry Smith PetscFunctionReturn(0); 3454b9ad928SBarry Smith } 3464b9ad928SBarry Smith 3474b9ad928SBarry Smith EXTERN_C_BEGIN 3484b9ad928SBarry Smith #undef __FUNCT__ 3494b9ad928SBarry Smith #define __FUNCT__ "PCRedundantGetPC_Redundant" 350dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCRedundantGetPC_Redundant(PC pc,PC *innerpc) 3514b9ad928SBarry Smith { 3525f06b7aaSBarry Smith PetscErrorCode ierr; 3534b9ad928SBarry Smith PC_Redundant *red = (PC_Redundant*)pc->data; 3545f06b7aaSBarry Smith MPI_Comm comm,subcomm; 3555f06b7aaSBarry Smith const char *prefix; 3564b9ad928SBarry Smith 3574b9ad928SBarry Smith PetscFunctionBegin; 3585f06b7aaSBarry Smith if (!red->psubcomm) { 3595f06b7aaSBarry Smith ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 360638faf0bSHong Zhang ierr = PetscSubcommCreate(comm,red->nsubcomm,PETSC_SUBCOMM_INTERLACED,&red->psubcomm);CHKERRQ(ierr); 3615f06b7aaSBarry Smith ierr = PetscLogObjectMemory(pc,sizeof(PetscSubcomm));CHKERRQ(ierr); 3625f06b7aaSBarry Smith 3635f06b7aaSBarry Smith /* create a new PC that processors in each subcomm have copy of */ 3645f06b7aaSBarry Smith subcomm = red->psubcomm->comm; 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 } 3765f06b7aaSBarry Smith 3775f06b7aaSBarry Smith ierr = KSPGetPC(red->ksp,innerpc);CHKERRQ(ierr); 3784b9ad928SBarry Smith PetscFunctionReturn(0); 3794b9ad928SBarry Smith } 3804b9ad928SBarry Smith EXTERN_C_END 3814b9ad928SBarry Smith 3824b9ad928SBarry Smith #undef __FUNCT__ 3834b9ad928SBarry Smith #define __FUNCT__ "PCRedundantGetPC" 3844b9ad928SBarry Smith /*@ 3854b9ad928SBarry Smith PCRedundantGetPC - Gets the sequential PC created by the redundant PC. 3864b9ad928SBarry Smith 3874b9ad928SBarry Smith Not Collective 3884b9ad928SBarry Smith 3894b9ad928SBarry Smith Input Parameter: 3904b9ad928SBarry Smith . pc - the preconditioner context 3914b9ad928SBarry Smith 3924b9ad928SBarry Smith Output Parameter: 3934b9ad928SBarry Smith . innerpc - the sequential PC 3944b9ad928SBarry Smith 3954b9ad928SBarry Smith Level: advanced 3964b9ad928SBarry Smith 3974b9ad928SBarry Smith .keywords: PC, redundant solve 3984b9ad928SBarry Smith @*/ 399dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCRedundantGetPC(PC pc,PC *innerpc) 4004b9ad928SBarry Smith { 401dfbe8321SBarry Smith PetscErrorCode ierr,(*f)(PC,PC*); 4024b9ad928SBarry Smith 4034b9ad928SBarry Smith PetscFunctionBegin; 4040700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 4054482741eSBarry Smith PetscValidPointer(innerpc,2); 4064b9ad928SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCRedundantGetPC_C",(void (**)(void))&f);CHKERRQ(ierr); 4074b9ad928SBarry Smith if (f) { 4084b9ad928SBarry Smith ierr = (*f)(pc,innerpc);CHKERRQ(ierr); 4094b9ad928SBarry Smith } 4104b9ad928SBarry Smith PetscFunctionReturn(0); 4114b9ad928SBarry Smith } 4124b9ad928SBarry Smith 4134b9ad928SBarry Smith EXTERN_C_BEGIN 4144b9ad928SBarry Smith #undef __FUNCT__ 4154b9ad928SBarry Smith #define __FUNCT__ "PCRedundantGetOperators_Redundant" 416dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCRedundantGetOperators_Redundant(PC pc,Mat *mat,Mat *pmat) 4174b9ad928SBarry Smith { 4184b9ad928SBarry Smith PC_Redundant *red = (PC_Redundant*)pc->data; 4194b9ad928SBarry Smith 4204b9ad928SBarry Smith PetscFunctionBegin; 421b3804887SHong Zhang if (mat) *mat = red->pmats; 422b3804887SHong Zhang if (pmat) *pmat = red->pmats; 4234b9ad928SBarry Smith PetscFunctionReturn(0); 4244b9ad928SBarry Smith } 4254b9ad928SBarry Smith EXTERN_C_END 4264b9ad928SBarry Smith 4274b9ad928SBarry Smith #undef __FUNCT__ 4284b9ad928SBarry Smith #define __FUNCT__ "PCRedundantGetOperators" 4294b9ad928SBarry Smith /*@ 4304b9ad928SBarry Smith PCRedundantGetOperators - gets the sequential matrix and preconditioner matrix 4314b9ad928SBarry Smith 4324b9ad928SBarry Smith Not Collective 4334b9ad928SBarry Smith 4344b9ad928SBarry Smith Input Parameter: 4354b9ad928SBarry Smith . pc - the preconditioner context 4364b9ad928SBarry Smith 4374b9ad928SBarry Smith Output Parameters: 4384b9ad928SBarry Smith + mat - the matrix 4394b9ad928SBarry Smith - pmat - the (possibly different) preconditioner matrix 4404b9ad928SBarry Smith 4414b9ad928SBarry Smith Level: advanced 4424b9ad928SBarry Smith 4434b9ad928SBarry Smith .keywords: PC, redundant solve 4444b9ad928SBarry Smith @*/ 445dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCRedundantGetOperators(PC pc,Mat *mat,Mat *pmat) 4464b9ad928SBarry Smith { 447dfbe8321SBarry Smith PetscErrorCode ierr,(*f)(PC,Mat*,Mat*); 4484b9ad928SBarry Smith 4494b9ad928SBarry Smith PetscFunctionBegin; 4500700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 4514482741eSBarry Smith if (mat) PetscValidPointer(mat,2); 4524482741eSBarry Smith if (pmat) PetscValidPointer(pmat,3); 4534b9ad928SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCRedundantGetOperators_C",(void (**)(void))&f);CHKERRQ(ierr); 4544b9ad928SBarry Smith if (f) { 4554b9ad928SBarry Smith ierr = (*f)(pc,mat,pmat);CHKERRQ(ierr); 4564b9ad928SBarry Smith } 4574b9ad928SBarry Smith PetscFunctionReturn(0); 4584b9ad928SBarry Smith } 4594b9ad928SBarry Smith 4604b9ad928SBarry Smith /* -------------------------------------------------------------------------------------*/ 46137a17b4dSBarry Smith /*MC 462a98ce0f4SHong Zhang PCREDUNDANT - Runs a preconditioner for the entire problem on subgroups of processors 46337a17b4dSBarry Smith 46437a17b4dSBarry Smith Options for the redundant preconditioners can be set with -redundant_pc_xxx 46537a17b4dSBarry Smith 46609391456SBarry Smith Options Database: 4679b21d695SBarry Smith . -pc_redundant_number <n> - number of redundant solves, for example if you are using 64 MPI processes and 4689b21d695SBarry Smith use an n of 4 there will be 4 parallel solves each on 16 = 64/4 processes. 46909391456SBarry Smith 47037a17b4dSBarry Smith Level: intermediate 47137a17b4dSBarry Smith 47237a17b4dSBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PCRedundantSetScatter(), 4739b21d695SBarry Smith PCRedundantGetPC(), PCRedundantGetOperators(), PCRedundantSetNumber() 47437a17b4dSBarry Smith M*/ 47537a17b4dSBarry Smith 4764b9ad928SBarry Smith EXTERN_C_BEGIN 4774b9ad928SBarry Smith #undef __FUNCT__ 4784b9ad928SBarry Smith #define __FUNCT__ "PCCreate_Redundant" 479dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_Redundant(PC pc) 4804b9ad928SBarry Smith { 481dfbe8321SBarry Smith PetscErrorCode ierr; 4824b9ad928SBarry Smith PC_Redundant *red; 48369db28dcSHong Zhang PetscMPIInt size; 4843f457be1SHong Zhang 4854b9ad928SBarry Smith PetscFunctionBegin; 48638f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_Redundant,&red);CHKERRQ(ierr); 4877adad957SLisandro Dalcin ierr = MPI_Comm_size(((PetscObject)pc)->comm,&size);CHKERRQ(ierr); 48869db28dcSHong Zhang red->nsubcomm = size; 4894b9ad928SBarry Smith red->useparallelmat = PETSC_TRUE; 4901fbd8f88SHong Zhang pc->data = (void*)red; 4914b9ad928SBarry Smith 4924b9ad928SBarry Smith pc->ops->apply = PCApply_Redundant; 4934b9ad928SBarry Smith pc->ops->applytranspose = 0; 4944b9ad928SBarry Smith pc->ops->setup = PCSetUp_Redundant; 4954b9ad928SBarry Smith pc->ops->destroy = PCDestroy_Redundant; 4964b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_Redundant; 4974b9ad928SBarry Smith pc->ops->view = PCView_Redundant; 4984b9ad928SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCRedundantSetScatter_C","PCRedundantSetScatter_Redundant", 4994b9ad928SBarry Smith PCRedundantSetScatter_Redundant);CHKERRQ(ierr); 50009a6bc64SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCRedundantSetNumber_C","PCRedundantSetNumber_Redundant", 50109a6bc64SHong Zhang PCRedundantSetNumber_Redundant);CHKERRQ(ierr); 5024b9ad928SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCRedundantGetPC_C","PCRedundantGetPC_Redundant", 5034b9ad928SBarry Smith PCRedundantGetPC_Redundant);CHKERRQ(ierr); 5044b9ad928SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCRedundantGetOperators_C","PCRedundantGetOperators_Redundant", 5054b9ad928SBarry Smith PCRedundantGetOperators_Redundant);CHKERRQ(ierr); 5064b9ad928SBarry Smith PetscFunctionReturn(0); 5074b9ad928SBarry Smith } 5084b9ad928SBarry Smith EXTERN_C_END 509