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 514b9ad928SBarry Smith #undef __FUNCT__ 524b9ad928SBarry Smith #define __FUNCT__ "PCSetUp_Redundant" 536849ba73SBarry Smith static PetscErrorCode PCSetUp_Redundant(PC pc) 544b9ad928SBarry Smith { 554b9ad928SBarry Smith PC_Redundant *red = (PC_Redundant*)pc->data; 56dfbe8321SBarry Smith PetscErrorCode ierr; 5745fc02eaSBarry Smith PetscInt mstart,mend,mlocal,m,mlocal_sub,rstart_sub,rend_sub,mloc_sub; 5813f74950SBarry Smith PetscMPIInt size; 594b9ad928SBarry Smith MatReuse reuse = MAT_INITIAL_MATRIX; 604b9ad928SBarry Smith MatStructure str = DIFFERENT_NONZERO_PATTERN; 61ce94432eSBarry Smith MPI_Comm comm,subcomm; 6223ce1328SBarry Smith Vec vec; 633f457be1SHong Zhang PetscMPIInt subsize,subrank; 641fbd8f88SHong Zhang const char *prefix; 65b862ddfaSBarry Smith const PetscInt *range; 663f457be1SHong Zhang 674b9ad928SBarry Smith PetscFunctionBegin; 68ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 6923ce1328SBarry Smith ierr = MatGetVecs(pc->pmat,&vec,0);CHKERRQ(ierr); 7023ce1328SBarry Smith ierr = VecGetSize(vec,&m);CHKERRQ(ierr); 711fbd8f88SHong Zhang 724b9ad928SBarry Smith if (!pc->setupcalled) { 735f06b7aaSBarry Smith if (!red->psubcomm) { 74d8a68f86SHong Zhang ierr = PetscSubcommCreate(comm,&red->psubcomm);CHKERRQ(ierr); 75d8a68f86SHong Zhang ierr = PetscSubcommSetNumber(red->psubcomm,red->nsubcomm);CHKERRQ(ierr); 76d8a68f86SHong Zhang ierr = PetscSubcommSetType(red->psubcomm,PETSC_SUBCOMM_INTERLACED);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; 812fa5cd67SKarl Rupp 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); 922fa5cd67SKarl Rupp } else subcomm = red->psubcomm->comm; 931fbd8f88SHong Zhang 943f457be1SHong Zhang /* create working vectors xsub/ysub and xdup/ydup */ 9523ce1328SBarry Smith ierr = VecGetLocalSize(vec,&mlocal);CHKERRQ(ierr); 963f457be1SHong Zhang ierr = VecGetOwnershipRange(vec,&mstart,&mend);CHKERRQ(ierr); 974b9ad928SBarry Smith 983f457be1SHong Zhang /* get local size of xsub/ysub */ 991fbd8f88SHong Zhang ierr = MPI_Comm_size(subcomm,&subsize);CHKERRQ(ierr); 1001fbd8f88SHong Zhang ierr = MPI_Comm_rank(subcomm,&subrank);CHKERRQ(ierr); 101b862ddfaSBarry Smith ierr = MatGetOwnershipRanges(pc->pmat,&range);CHKERRQ(ierr); 102b862ddfaSBarry Smith rstart_sub = range[red->psubcomm->n*subrank]; /* rstart in xsub/ysub */ 1032fa5cd67SKarl Rupp if (subrank+1 < subsize) rend_sub = range[red->psubcomm->n*(subrank+1)]; 1042fa5cd67SKarl Rupp else rend_sub = m; 1052fa5cd67SKarl Rupp 1063f457be1SHong Zhang mloc_sub = rend_sub - rstart_sub; 1071fbd8f88SHong Zhang ierr = VecCreateMPI(subcomm,mloc_sub,PETSC_DECIDE,&red->ysub);CHKERRQ(ierr); 1083f457be1SHong Zhang /* create xsub with empty local arrays, because xdup's arrays will be placed into it */ 1090298fd71SBarry Smith ierr = VecCreateMPIWithArray(subcomm,1,mloc_sub,PETSC_DECIDE,NULL,&red->xsub);CHKERRQ(ierr); 1103f457be1SHong Zhang 1113f457be1SHong Zhang /* create xdup and ydup. ydup has empty local arrays because ysub's arrays will be place into it. 112ce94432eSBarry Smith Note: we use communicator dupcomm, not PetscObjectComm((PetscObject)pc)! */ 1131fbd8f88SHong Zhang ierr = VecCreateMPI(red->psubcomm->dupparent,mloc_sub,PETSC_DECIDE,&red->xdup);CHKERRQ(ierr); 1140298fd71SBarry Smith ierr = VecCreateMPIWithArray(red->psubcomm->dupparent,1,mloc_sub,PETSC_DECIDE,NULL,&red->ydup);CHKERRQ(ierr); 1153f457be1SHong Zhang 1163f457be1SHong Zhang /* create vec scatters */ 1173f457be1SHong Zhang if (!red->scatterin) { 1183f457be1SHong Zhang IS is1,is2; 1193f457be1SHong Zhang PetscInt *idx1,*idx2,i,j,k; 12045fc02eaSBarry Smith 1211d79065fSBarry Smith ierr = PetscMalloc2(red->psubcomm->n*mlocal,PetscInt,&idx1,red->psubcomm->n*mlocal,PetscInt,&idx2);CHKERRQ(ierr); 1223f457be1SHong Zhang j = 0; 1231fbd8f88SHong Zhang for (k=0; k<red->psubcomm->n; k++) { 1243f457be1SHong Zhang for (i=mstart; i<mend; i++) { 1253f457be1SHong Zhang idx1[j] = i; 1263f457be1SHong Zhang idx2[j++] = i + m*k; 1273f457be1SHong Zhang } 1283f457be1SHong Zhang } 12970b3c8c7SBarry Smith ierr = ISCreateGeneral(comm,red->psubcomm->n*mlocal,idx1,PETSC_COPY_VALUES,&is1);CHKERRQ(ierr); 13070b3c8c7SBarry Smith ierr = ISCreateGeneral(comm,red->psubcomm->n*mlocal,idx2,PETSC_COPY_VALUES,&is2);CHKERRQ(ierr); 1313f457be1SHong Zhang ierr = VecScatterCreate(vec,is1,red->xdup,is2,&red->scatterin);CHKERRQ(ierr); 132fcfd50ebSBarry Smith ierr = ISDestroy(&is1);CHKERRQ(ierr); 133fcfd50ebSBarry Smith ierr = ISDestroy(&is2);CHKERRQ(ierr); 1343f457be1SHong Zhang 1351fbd8f88SHong Zhang ierr = ISCreateStride(comm,mlocal,mstart+ red->psubcomm->color*m,1,&is1);CHKERRQ(ierr); 1363f457be1SHong Zhang ierr = ISCreateStride(comm,mlocal,mstart,1,&is2);CHKERRQ(ierr); 1373f457be1SHong Zhang ierr = VecScatterCreate(red->xdup,is1,vec,is2,&red->scatterout);CHKERRQ(ierr); 138fcfd50ebSBarry Smith ierr = ISDestroy(&is1);CHKERRQ(ierr); 139fcfd50ebSBarry Smith ierr = ISDestroy(&is2);CHKERRQ(ierr); 1401d79065fSBarry Smith ierr = PetscFree2(idx1,idx2);CHKERRQ(ierr); 1414b9ad928SBarry Smith } 1424b9ad928SBarry Smith } 1436bf464f9SBarry Smith ierr = VecDestroy(&vec);CHKERRQ(ierr); 1444b9ad928SBarry Smith 1454b9ad928SBarry Smith /* if pmatrix set by user is sequential then we do not need to gather the parallel matrix */ 1463f457be1SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1472fa5cd67SKarl Rupp if (size == 1) red->useparallelmat = PETSC_FALSE; 1484b9ad928SBarry Smith 1494b9ad928SBarry Smith if (red->useparallelmat) { 1504b9ad928SBarry Smith if (pc->setupcalled == 1 && pc->flag == DIFFERENT_NONZERO_PATTERN) { 1514b9ad928SBarry Smith /* destroy old matrices */ 1526bf464f9SBarry Smith ierr = MatDestroy(&red->pmats);CHKERRQ(ierr); 1534b9ad928SBarry Smith } else if (pc->setupcalled == 1) { 1544b9ad928SBarry Smith reuse = MAT_REUSE_MATRIX; 1554b9ad928SBarry Smith str = SAME_NONZERO_PATTERN; 1564b9ad928SBarry Smith } 1574b9ad928SBarry Smith 1583f457be1SHong Zhang /* grab the parallel matrix and put it into processors of a subcomminicator */ 159f664ae05SHong Zhang /*--------------------------------------------------------------------------*/ 160f664ae05SHong Zhang ierr = VecGetLocalSize(red->ysub,&mlocal_sub);CHKERRQ(ierr); 16169db28dcSHong Zhang ierr = MatGetRedundantMatrix(pc->pmat,red->psubcomm->n,red->psubcomm->comm,mlocal_sub,reuse,&red->pmats);CHKERRQ(ierr); 1623f457be1SHong Zhang /* tell PC of the subcommunicator its operators */ 16390f1c854SHong Zhang ierr = KSPSetOperators(red->ksp,red->pmats,red->pmats,str);CHKERRQ(ierr); 1644b9ad928SBarry Smith } else { 16590f1c854SHong Zhang ierr = KSPSetOperators(red->ksp,pc->mat,pc->pmat,pc->flag);CHKERRQ(ierr); 1664b9ad928SBarry Smith } 1670c24e6a1SHong Zhang if (pc->setfromoptionscalled) { 1683e065800SHong Zhang ierr = KSPSetFromOptions(red->ksp);CHKERRQ(ierr); 1690c24e6a1SHong Zhang } 1703e065800SHong Zhang ierr = KSPSetUp(red->ksp);CHKERRQ(ierr); 1714b9ad928SBarry Smith PetscFunctionReturn(0); 1724b9ad928SBarry Smith } 1734b9ad928SBarry Smith 1744b9ad928SBarry Smith #undef __FUNCT__ 1754b9ad928SBarry Smith #define __FUNCT__ "PCApply_Redundant" 1766849ba73SBarry Smith static PetscErrorCode PCApply_Redundant(PC pc,Vec x,Vec y) 1774b9ad928SBarry Smith { 1784b9ad928SBarry Smith PC_Redundant *red = (PC_Redundant*)pc->data; 179dfbe8321SBarry Smith PetscErrorCode ierr; 1803f457be1SHong Zhang PetscScalar *array; 1814b9ad928SBarry Smith 1824b9ad928SBarry Smith PetscFunctionBegin; 1833f457be1SHong Zhang /* scatter x to xdup */ 184ca9f406cSSatish Balay ierr = VecScatterBegin(red->scatterin,x,red->xdup,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 185ca9f406cSSatish Balay ierr = VecScatterEnd(red->scatterin,x,red->xdup,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1863f457be1SHong Zhang 1873f457be1SHong Zhang /* place xdup's local array into xsub */ 1883f457be1SHong Zhang ierr = VecGetArray(red->xdup,&array);CHKERRQ(ierr); 1893f457be1SHong Zhang ierr = VecPlaceArray(red->xsub,(const PetscScalar*)array);CHKERRQ(ierr); 1904b9ad928SBarry Smith 1914b9ad928SBarry Smith /* apply preconditioner on each processor */ 19283ab6a24SBarry Smith ierr = KSPSolve(red->ksp,red->xsub,red->ysub);CHKERRQ(ierr); 1933f457be1SHong Zhang ierr = VecResetArray(red->xsub);CHKERRQ(ierr); 1943f457be1SHong Zhang ierr = VecRestoreArray(red->xdup,&array);CHKERRQ(ierr); 1954b9ad928SBarry Smith 1963f457be1SHong Zhang /* place ysub's local array into ydup */ 1973f457be1SHong Zhang ierr = VecGetArray(red->ysub,&array);CHKERRQ(ierr); 1983f457be1SHong Zhang ierr = VecPlaceArray(red->ydup,(const PetscScalar*)array);CHKERRQ(ierr); 1993f457be1SHong Zhang 2003f457be1SHong Zhang /* scatter ydup to y */ 201ca9f406cSSatish Balay ierr = VecScatterBegin(red->scatterout,red->ydup,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 202ca9f406cSSatish Balay ierr = VecScatterEnd(red->scatterout,red->ydup,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 2033f457be1SHong Zhang ierr = VecResetArray(red->ydup);CHKERRQ(ierr); 2043f457be1SHong Zhang ierr = VecRestoreArray(red->ysub,&array);CHKERRQ(ierr); 2054b9ad928SBarry Smith PetscFunctionReturn(0); 2064b9ad928SBarry Smith } 2074b9ad928SBarry Smith 2084b9ad928SBarry Smith #undef __FUNCT__ 2091ea5a559SBarry Smith #define __FUNCT__ "PCReset_Redundant" 2101ea5a559SBarry Smith static PetscErrorCode PCReset_Redundant(PC pc) 2114b9ad928SBarry Smith { 2124b9ad928SBarry Smith PC_Redundant *red = (PC_Redundant*)pc->data; 213dfbe8321SBarry Smith PetscErrorCode ierr; 2144b9ad928SBarry Smith 2154b9ad928SBarry Smith PetscFunctionBegin; 2166bf464f9SBarry Smith ierr = VecScatterDestroy(&red->scatterin);CHKERRQ(ierr); 2176bf464f9SBarry Smith ierr = VecScatterDestroy(&red->scatterout);CHKERRQ(ierr); 2186bf464f9SBarry Smith ierr = VecDestroy(&red->ysub);CHKERRQ(ierr); 2196bf464f9SBarry Smith ierr = VecDestroy(&red->xsub);CHKERRQ(ierr); 2206bf464f9SBarry Smith ierr = VecDestroy(&red->xdup);CHKERRQ(ierr); 2216bf464f9SBarry Smith ierr = VecDestroy(&red->ydup);CHKERRQ(ierr); 2226bf464f9SBarry Smith ierr = MatDestroy(&red->pmats);CHKERRQ(ierr); 223546c30cbSLisandro Dalcin if (red->ksp) {ierr = KSPReset(red->ksp);CHKERRQ(ierr);} 2241ea5a559SBarry Smith PetscFunctionReturn(0); 2251ea5a559SBarry Smith } 2261ea5a559SBarry Smith 2271ea5a559SBarry Smith #undef __FUNCT__ 2281ea5a559SBarry Smith #define __FUNCT__ "PCDestroy_Redundant" 2291ea5a559SBarry Smith static PetscErrorCode PCDestroy_Redundant(PC pc) 2301ea5a559SBarry Smith { 2311ea5a559SBarry Smith PC_Redundant *red = (PC_Redundant*)pc->data; 2321ea5a559SBarry Smith PetscErrorCode ierr; 2331ea5a559SBarry Smith 2341ea5a559SBarry Smith PetscFunctionBegin; 2351ea5a559SBarry Smith ierr = PCReset_Redundant(pc);CHKERRQ(ierr); 2366bf464f9SBarry Smith ierr = KSPDestroy(&red->ksp);CHKERRQ(ierr); 2376bf464f9SBarry Smith ierr = PetscSubcommDestroy(&red->psubcomm);CHKERRQ(ierr); 238c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 2394b9ad928SBarry Smith PetscFunctionReturn(0); 2404b9ad928SBarry Smith } 2414b9ad928SBarry Smith 2424b9ad928SBarry Smith #undef __FUNCT__ 2434b9ad928SBarry Smith #define __FUNCT__ "PCSetFromOptions_Redundant" 2446849ba73SBarry Smith static PetscErrorCode PCSetFromOptions_Redundant(PC pc) 2454b9ad928SBarry Smith { 246a98ce0f4SHong Zhang PetscErrorCode ierr; 247a98ce0f4SHong Zhang PC_Redundant *red = (PC_Redundant*)pc->data; 248a98ce0f4SHong Zhang 2494b9ad928SBarry Smith PetscFunctionBegin; 250a98ce0f4SHong Zhang ierr = PetscOptionsHead("Redundant options");CHKERRQ(ierr); 25109a6bc64SHong Zhang ierr = PetscOptionsInt("-pc_redundant_number","Number of redundant pc","PCRedundantSetNumber",red->nsubcomm,&red->nsubcomm,0);CHKERRQ(ierr); 252a98ce0f4SHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 2534b9ad928SBarry Smith PetscFunctionReturn(0); 2544b9ad928SBarry Smith } 2554b9ad928SBarry Smith 2564b9ad928SBarry Smith #undef __FUNCT__ 25709a6bc64SHong Zhang #define __FUNCT__ "PCRedundantSetNumber_Redundant" 258f7a08781SBarry Smith static PetscErrorCode PCRedundantSetNumber_Redundant(PC pc,PetscInt nreds) 25909a6bc64SHong Zhang { 26009a6bc64SHong Zhang PC_Redundant *red = (PC_Redundant*)pc->data; 26109a6bc64SHong Zhang 26209a6bc64SHong Zhang PetscFunctionBegin; 26309a6bc64SHong Zhang red->nsubcomm = nreds; 26409a6bc64SHong Zhang PetscFunctionReturn(0); 26509a6bc64SHong Zhang } 26609a6bc64SHong Zhang 26709a6bc64SHong Zhang #undef __FUNCT__ 26809a6bc64SHong Zhang #define __FUNCT__ "PCRedundantSetNumber" 26909a6bc64SHong Zhang /*@ 27009a6bc64SHong Zhang PCRedundantSetNumber - Sets the number of redundant preconditioner contexts. 27109a6bc64SHong Zhang 2723f9fe445SBarry Smith Logically Collective on PC 27309a6bc64SHong Zhang 27409a6bc64SHong Zhang Input Parameters: 27509a6bc64SHong Zhang + pc - the preconditioner context 2769b21d695SBarry Smith - nredundant - number of redundant preconditioner contexts; for example if you are using 64 MPI processes and 2779b21d695SBarry Smith use an nredundant of 4 there will be 4 parallel solves each on 16 = 64/4 processes. 27809a6bc64SHong Zhang 27909a6bc64SHong Zhang Level: advanced 28009a6bc64SHong Zhang 28109a6bc64SHong Zhang .keywords: PC, redundant solve 28209a6bc64SHong Zhang @*/ 2837087cfbeSBarry Smith PetscErrorCode PCRedundantSetNumber(PC pc,PetscInt nredundant) 28409a6bc64SHong Zhang { 2854ac538c5SBarry Smith PetscErrorCode ierr; 28609a6bc64SHong Zhang 28709a6bc64SHong Zhang PetscFunctionBegin; 2880700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 289ce94432eSBarry Smith if (nredundant <= 0) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG, "num of redundant pc %D must be positive",nredundant); 2904ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCRedundantSetNumber_C",(PC,PetscInt),(pc,nredundant));CHKERRQ(ierr); 29109a6bc64SHong Zhang PetscFunctionReturn(0); 29209a6bc64SHong Zhang } 29309a6bc64SHong Zhang 29409a6bc64SHong Zhang #undef __FUNCT__ 2954b9ad928SBarry Smith #define __FUNCT__ "PCRedundantSetScatter_Redundant" 296f7a08781SBarry Smith static PetscErrorCode PCRedundantSetScatter_Redundant(PC pc,VecScatter in,VecScatter out) 2974b9ad928SBarry Smith { 2984b9ad928SBarry Smith PC_Redundant *red = (PC_Redundant*)pc->data; 299dfbe8321SBarry Smith PetscErrorCode ierr; 3004b9ad928SBarry Smith 3014b9ad928SBarry Smith PetscFunctionBegin; 3024b9ad928SBarry Smith ierr = PetscObjectReference((PetscObject)in);CHKERRQ(ierr); 3036bf464f9SBarry Smith ierr = VecScatterDestroy(&red->scatterin);CHKERRQ(ierr); 3042fa5cd67SKarl Rupp 305c3122656SLisandro Dalcin red->scatterin = in; 3062fa5cd67SKarl Rupp 3074b9ad928SBarry Smith ierr = PetscObjectReference((PetscObject)out);CHKERRQ(ierr); 3086bf464f9SBarry Smith ierr = VecScatterDestroy(&red->scatterout);CHKERRQ(ierr); 309c3122656SLisandro Dalcin red->scatterout = out; 3104b9ad928SBarry Smith PetscFunctionReturn(0); 3114b9ad928SBarry Smith } 3124b9ad928SBarry Smith 3134b9ad928SBarry Smith #undef __FUNCT__ 3144b9ad928SBarry Smith #define __FUNCT__ "PCRedundantSetScatter" 3154b9ad928SBarry Smith /*@ 3164b9ad928SBarry Smith PCRedundantSetScatter - Sets the scatter used to copy values into the 3174b9ad928SBarry Smith redundant local solve and the scatter to move them back into the global 3184b9ad928SBarry Smith vector. 3194b9ad928SBarry Smith 3203f9fe445SBarry Smith Logically Collective on PC 3214b9ad928SBarry Smith 3224b9ad928SBarry Smith Input Parameters: 3234b9ad928SBarry Smith + pc - the preconditioner context 3244b9ad928SBarry Smith . in - the scatter to move the values in 3254b9ad928SBarry Smith - out - the scatter to move them out 3264b9ad928SBarry Smith 3274b9ad928SBarry Smith Level: advanced 3284b9ad928SBarry Smith 3294b9ad928SBarry Smith .keywords: PC, redundant solve 3304b9ad928SBarry Smith @*/ 3317087cfbeSBarry Smith PetscErrorCode PCRedundantSetScatter(PC pc,VecScatter in,VecScatter out) 3324b9ad928SBarry Smith { 3334ac538c5SBarry Smith PetscErrorCode ierr; 3344b9ad928SBarry Smith 3354b9ad928SBarry Smith PetscFunctionBegin; 3360700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 3370700a824SBarry Smith PetscValidHeaderSpecific(in,VEC_SCATTER_CLASSID,2); 3380700a824SBarry Smith PetscValidHeaderSpecific(out,VEC_SCATTER_CLASSID,3); 3394ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCRedundantSetScatter_C",(PC,VecScatter,VecScatter),(pc,in,out));CHKERRQ(ierr); 3404b9ad928SBarry Smith PetscFunctionReturn(0); 3414b9ad928SBarry Smith } 3424b9ad928SBarry Smith 3434b9ad928SBarry Smith #undef __FUNCT__ 34483ab6a24SBarry Smith #define __FUNCT__ "PCRedundantGetKSP_Redundant" 345f7a08781SBarry Smith static PetscErrorCode PCRedundantGetKSP_Redundant(PC pc,KSP *innerksp) 3464b9ad928SBarry Smith { 3475f06b7aaSBarry Smith PetscErrorCode ierr; 3484b9ad928SBarry Smith PC_Redundant *red = (PC_Redundant*)pc->data; 3495f06b7aaSBarry Smith MPI_Comm comm,subcomm; 3505f06b7aaSBarry Smith const char *prefix; 3514b9ad928SBarry Smith 3524b9ad928SBarry Smith PetscFunctionBegin; 3535f06b7aaSBarry Smith if (!red->psubcomm) { 3545f06b7aaSBarry Smith ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 355d8a68f86SHong Zhang ierr = PetscSubcommCreate(comm,&red->psubcomm);CHKERRQ(ierr); 356d8a68f86SHong Zhang ierr = PetscSubcommSetNumber(red->psubcomm,red->nsubcomm);CHKERRQ(ierr); 357d8a68f86SHong Zhang ierr = PetscSubcommSetType(red->psubcomm,PETSC_SUBCOMM_INTERLACED);CHKERRQ(ierr); 3585f06b7aaSBarry Smith ierr = PetscLogObjectMemory(pc,sizeof(PetscSubcomm));CHKERRQ(ierr); 3595f06b7aaSBarry Smith 3605f06b7aaSBarry Smith /* create a new PC that processors in each subcomm have copy of */ 3615f06b7aaSBarry Smith subcomm = red->psubcomm->comm; 3622fa5cd67SKarl Rupp 3635f06b7aaSBarry Smith ierr = KSPCreate(subcomm,&red->ksp);CHKERRQ(ierr); 3645f06b7aaSBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)red->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 3655f06b7aaSBarry Smith ierr = PetscLogObjectParent(pc,red->ksp);CHKERRQ(ierr); 3665f06b7aaSBarry Smith ierr = KSPSetType(red->ksp,KSPPREONLY);CHKERRQ(ierr); 3675f06b7aaSBarry Smith ierr = KSPGetPC(red->ksp,&red->pc);CHKERRQ(ierr); 3685f06b7aaSBarry Smith ierr = PCSetType(red->pc,PCLU);CHKERRQ(ierr); 3695f06b7aaSBarry Smith 3705f06b7aaSBarry Smith ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 3715f06b7aaSBarry Smith ierr = KSPSetOptionsPrefix(red->ksp,prefix);CHKERRQ(ierr); 3725f06b7aaSBarry Smith ierr = KSPAppendOptionsPrefix(red->ksp,"redundant_");CHKERRQ(ierr); 3735f06b7aaSBarry Smith } 37483ab6a24SBarry Smith *innerksp = red->ksp; 3754b9ad928SBarry Smith PetscFunctionReturn(0); 3764b9ad928SBarry Smith } 3774b9ad928SBarry Smith 3784b9ad928SBarry Smith #undef __FUNCT__ 37983ab6a24SBarry Smith #define __FUNCT__ "PCRedundantGetKSP" 3804b9ad928SBarry Smith /*@ 38183ab6a24SBarry Smith PCRedundantGetKSP - Gets the less parallel KSP created by the redundant PC. 3824b9ad928SBarry Smith 3834b9ad928SBarry Smith Not Collective 3844b9ad928SBarry Smith 3854b9ad928SBarry Smith Input Parameter: 3864b9ad928SBarry Smith . pc - the preconditioner context 3874b9ad928SBarry Smith 3884b9ad928SBarry Smith Output Parameter: 38983ab6a24SBarry Smith . innerksp - the KSP on the smaller set of processes 3904b9ad928SBarry Smith 3914b9ad928SBarry Smith Level: advanced 3924b9ad928SBarry Smith 3934b9ad928SBarry Smith .keywords: PC, redundant solve 3944b9ad928SBarry Smith @*/ 39583ab6a24SBarry Smith PetscErrorCode PCRedundantGetKSP(PC pc,KSP *innerksp) 3964b9ad928SBarry Smith { 3974ac538c5SBarry Smith PetscErrorCode ierr; 3984b9ad928SBarry Smith 3994b9ad928SBarry Smith PetscFunctionBegin; 4000700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 40183ab6a24SBarry Smith PetscValidPointer(innerksp,2); 40283ab6a24SBarry Smith ierr = PetscTryMethod(pc,"PCRedundantGetKSP_C",(PC,KSP*),(pc,innerksp));CHKERRQ(ierr); 4034b9ad928SBarry Smith PetscFunctionReturn(0); 4044b9ad928SBarry Smith } 4054b9ad928SBarry Smith 4064b9ad928SBarry Smith #undef __FUNCT__ 4074b9ad928SBarry Smith #define __FUNCT__ "PCRedundantGetOperators_Redundant" 408f7a08781SBarry Smith static PetscErrorCode PCRedundantGetOperators_Redundant(PC pc,Mat *mat,Mat *pmat) 4094b9ad928SBarry Smith { 4104b9ad928SBarry Smith PC_Redundant *red = (PC_Redundant*)pc->data; 4114b9ad928SBarry Smith 4124b9ad928SBarry Smith PetscFunctionBegin; 413b3804887SHong Zhang if (mat) *mat = red->pmats; 414b3804887SHong Zhang if (pmat) *pmat = red->pmats; 4154b9ad928SBarry Smith PetscFunctionReturn(0); 4164b9ad928SBarry Smith } 4174b9ad928SBarry Smith 4184b9ad928SBarry Smith #undef __FUNCT__ 4194b9ad928SBarry Smith #define __FUNCT__ "PCRedundantGetOperators" 4204b9ad928SBarry Smith /*@ 4214b9ad928SBarry Smith PCRedundantGetOperators - gets the sequential matrix and preconditioner matrix 4224b9ad928SBarry Smith 4234b9ad928SBarry Smith Not Collective 4244b9ad928SBarry Smith 4254b9ad928SBarry Smith Input Parameter: 4264b9ad928SBarry Smith . pc - the preconditioner context 4274b9ad928SBarry Smith 4284b9ad928SBarry Smith Output Parameters: 4294b9ad928SBarry Smith + mat - the matrix 4304b9ad928SBarry Smith - pmat - the (possibly different) preconditioner matrix 4314b9ad928SBarry Smith 4324b9ad928SBarry Smith Level: advanced 4334b9ad928SBarry Smith 4344b9ad928SBarry Smith .keywords: PC, redundant solve 4354b9ad928SBarry Smith @*/ 4367087cfbeSBarry Smith PetscErrorCode PCRedundantGetOperators(PC pc,Mat *mat,Mat *pmat) 4374b9ad928SBarry Smith { 4384ac538c5SBarry Smith PetscErrorCode ierr; 4394b9ad928SBarry Smith 4404b9ad928SBarry Smith PetscFunctionBegin; 4410700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 4424482741eSBarry Smith if (mat) PetscValidPointer(mat,2); 4434482741eSBarry Smith if (pmat) PetscValidPointer(pmat,3); 4444ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCRedundantGetOperators_C",(PC,Mat*,Mat*),(pc,mat,pmat));CHKERRQ(ierr); 4454b9ad928SBarry Smith PetscFunctionReturn(0); 4464b9ad928SBarry Smith } 4474b9ad928SBarry Smith 4484b9ad928SBarry Smith /* -------------------------------------------------------------------------------------*/ 44937a17b4dSBarry Smith /*MC 45083ab6a24SBarry Smith PCREDUNDANT - Runs a KSP solver with preconditioner for the entire problem on subgroups of processors 45137a17b4dSBarry Smith 45283ab6a24SBarry Smith Options for the redundant preconditioners can be set with -redundant_pc_xxx for the redundant KSP with -redundant_ksp_xxx 45337a17b4dSBarry Smith 45409391456SBarry Smith Options Database: 4559b21d695SBarry Smith . -pc_redundant_number <n> - number of redundant solves, for example if you are using 64 MPI processes and 4569b21d695SBarry Smith use an n of 4 there will be 4 parallel solves each on 16 = 64/4 processes. 45709391456SBarry Smith 45837a17b4dSBarry Smith Level: intermediate 45937a17b4dSBarry Smith 46083ab6a24SBarry Smith Notes: The default KSP is preonly and the default PC is LU. 46183ab6a24SBarry Smith 46283ab6a24SBarry Smith Developer Notes: Note that PCSetInitialGuessNonzero() is not used by this class but likely should be. 4639cfaa89bSBarry Smith 46437a17b4dSBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PCRedundantSetScatter(), 46583ab6a24SBarry Smith PCRedundantGetKSP(), PCRedundantGetOperators(), PCRedundantSetNumber() 46637a17b4dSBarry Smith M*/ 46737a17b4dSBarry Smith 4684b9ad928SBarry Smith #undef __FUNCT__ 4694b9ad928SBarry Smith #define __FUNCT__ "PCCreate_Redundant" 470*8cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_Redundant(PC pc) 4714b9ad928SBarry Smith { 472dfbe8321SBarry Smith PetscErrorCode ierr; 4734b9ad928SBarry Smith PC_Redundant *red; 47469db28dcSHong Zhang PetscMPIInt size; 4753f457be1SHong Zhang 4764b9ad928SBarry Smith PetscFunctionBegin; 47738f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_Redundant,&red);CHKERRQ(ierr); 478ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr); 4792fa5cd67SKarl Rupp 48069db28dcSHong Zhang red->nsubcomm = size; 4814b9ad928SBarry Smith red->useparallelmat = PETSC_TRUE; 4821fbd8f88SHong Zhang pc->data = (void*)red; 4834b9ad928SBarry Smith 4844b9ad928SBarry Smith pc->ops->apply = PCApply_Redundant; 4854b9ad928SBarry Smith pc->ops->applytranspose = 0; 4864b9ad928SBarry Smith pc->ops->setup = PCSetUp_Redundant; 4874b9ad928SBarry Smith pc->ops->destroy = PCDestroy_Redundant; 4881ea5a559SBarry Smith pc->ops->reset = PCReset_Redundant; 4894b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_Redundant; 4904b9ad928SBarry Smith pc->ops->view = PCView_Redundant; 4912fa5cd67SKarl Rupp 49200de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCRedundantSetScatter_C","PCRedundantSetScatter_Redundant",PCRedundantSetScatter_Redundant);CHKERRQ(ierr); 49300de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCRedundantSetNumber_C","PCRedundantSetNumber_Redundant",PCRedundantSetNumber_Redundant);CHKERRQ(ierr); 49400de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCRedundantGetKSP_C","PCRedundantGetKSP_Redundant",PCRedundantGetKSP_Redundant);CHKERRQ(ierr); 49500de8ff0SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCRedundantGetOperators_C","PCRedundantGetOperators_Redundant",PCRedundantGetOperators_Redundant);CHKERRQ(ierr); 4964b9ad928SBarry Smith PetscFunctionReturn(0); 4974b9ad928SBarry Smith } 498b2573a8aSBarry Smith 499