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 51d3b23db5SHong Zhang extern PetscErrorCode MatGetRedundantMatrix_MPIAIJ_interlaced(Mat,PetscInt,PetscSubcomm,MatReuse,Mat*); /* rm later! */ 52d3b23db5SHong 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; 591b81debcSHong Zhang PetscInt mstart,mend,mlocal,M; 6013f74950SBarry Smith PetscMPIInt size; 61ce94432eSBarry Smith MPI_Comm comm,subcomm; 62ddc54837SHong Zhang Vec x; 631fbd8f88SHong Zhang const char *prefix; 643f457be1SHong Zhang 654b9ad928SBarry Smith PetscFunctionBegin; 66ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 67ddc54837SHong Zhang 68ddc54837SHong Zhang /* if pmatrix set by user is sequential then we do not need to gather the parallel matrix */ 69ddc54837SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 70ddc54837SHong Zhang if (size == 1) red->useparallelmat = PETSC_FALSE; 711fbd8f88SHong Zhang 724b9ad928SBarry Smith if (!pc->setupcalled) { 731b81debcSHong Zhang PetscInt mloc_sub; 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; 822fa5cd67SKarl Rupp 835f06b7aaSBarry Smith ierr = KSPCreate(subcomm,&red->ksp);CHKERRQ(ierr); 845f06b7aaSBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)red->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 855f06b7aaSBarry Smith ierr = PetscLogObjectParent(pc,red->ksp);CHKERRQ(ierr); 865f06b7aaSBarry Smith ierr = KSPSetType(red->ksp,KSPPREONLY);CHKERRQ(ierr); 875f06b7aaSBarry Smith ierr = KSPGetPC(red->ksp,&red->pc);CHKERRQ(ierr); 88cf52b8b1SHong Zhang ierr = PCSetType(red->pc,PCLU);CHKERRQ(ierr); 89cf52b8b1SHong Zhang 901fbd8f88SHong Zhang ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 915f06b7aaSBarry Smith ierr = KSPSetOptionsPrefix(red->ksp,prefix);CHKERRQ(ierr); 925f06b7aaSBarry Smith ierr = KSPAppendOptionsPrefix(red->ksp,"redundant_");CHKERRQ(ierr); 931b81debcSHong Zhang } else { 941b81debcSHong Zhang subcomm = red->psubcomm->comm; 951b81debcSHong Zhang } 961fbd8f88SHong Zhang 971b81debcSHong Zhang if (red->useparallelmat) { 981b81debcSHong Zhang /* grab the parallel matrix and put it into processors of a subcomminicator */ 991b81debcSHong Zhang ierr = MatGetRedundantMatrix_MPIAIJ_interlaced(pc->pmat,red->psubcomm->n,red->psubcomm,MAT_INITIAL_MATRIX,&red->pmats);CHKERRQ(ierr); 1001b81debcSHong Zhang ierr = KSPSetOperators(red->ksp,red->pmats,red->pmats,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 1014b9ad928SBarry Smith 1021b81debcSHong Zhang /* get working vectors xsub and ysub */ 1031b81debcSHong Zhang ierr = MatGetVecs(red->pmats,&red->xsub,&red->ysub);CHKERRQ(ierr); 1042fa5cd67SKarl Rupp 1051b81debcSHong Zhang /* create working vectors xdup and ydup. ydup has empty local arrays because ysub's arrays will be place into it. 106ce94432eSBarry Smith Note: we use communicator dupcomm, not PetscObjectComm((PetscObject)pc)! */ 1071b81debcSHong Zhang ierr = MatGetLocalSize(red->pmats,&mloc_sub,NULL);CHKERRQ(ierr); 1081fbd8f88SHong Zhang ierr = VecCreateMPI(red->psubcomm->dupparent,mloc_sub,PETSC_DECIDE,&red->xdup);CHKERRQ(ierr); 1090298fd71SBarry Smith ierr = VecCreateMPIWithArray(red->psubcomm->dupparent,1,mloc_sub,PETSC_DECIDE,NULL,&red->ydup);CHKERRQ(ierr); 1103f457be1SHong Zhang 111ab661555SHong Zhang /* create vecscatters for PETSC_SUBCOMM_INTERLACED! */ 1123f457be1SHong Zhang if (!red->scatterin) { 1133f457be1SHong Zhang IS is1,is2; 1143f457be1SHong Zhang PetscInt *idx1,*idx2,i,j,k; 11545fc02eaSBarry Smith 1161b81debcSHong Zhang ierr = MatGetVecs(pc->pmat,&x,0);CHKERRQ(ierr); 1171b81debcSHong Zhang ierr = VecGetSize(x,&M);CHKERRQ(ierr); 1181b81debcSHong Zhang ierr = VecGetOwnershipRange(x,&mstart,&mend);CHKERRQ(ierr); 1191b81debcSHong Zhang mlocal = mend - mstart; 120*1714dc9eSHong Zhang ierr = PetscMalloc2(red->psubcomm->n*mlocal,PetscInt,&idx1,red->psubcomm->n*mlocal,PetscInt,&idx2);CHKERRQ(ierr); 1213f457be1SHong Zhang j = 0; 1221fbd8f88SHong Zhang for (k=0; k<red->psubcomm->n; k++) { 1233f457be1SHong Zhang for (i=mstart; i<mend; i++) { 1243f457be1SHong Zhang idx1[j] = i; 125ddc54837SHong Zhang idx2[j++] = i + M*k; 1263f457be1SHong Zhang } 1273f457be1SHong Zhang } 12870b3c8c7SBarry Smith ierr = ISCreateGeneral(comm,red->psubcomm->n*mlocal,idx1,PETSC_COPY_VALUES,&is1);CHKERRQ(ierr); 12970b3c8c7SBarry Smith ierr = ISCreateGeneral(comm,red->psubcomm->n*mlocal,idx2,PETSC_COPY_VALUES,&is2);CHKERRQ(ierr); 130ddc54837SHong Zhang ierr = VecScatterCreate(x,is1,red->xdup,is2,&red->scatterin);CHKERRQ(ierr); 131fcfd50ebSBarry Smith ierr = ISDestroy(&is1);CHKERRQ(ierr); 132fcfd50ebSBarry Smith ierr = ISDestroy(&is2);CHKERRQ(ierr); 1333f457be1SHong Zhang 134ddc54837SHong Zhang ierr = ISCreateStride(comm,mlocal,mstart+ red->psubcomm->color*M,1,&is1);CHKERRQ(ierr); 1353f457be1SHong Zhang ierr = ISCreateStride(comm,mlocal,mstart,1,&is2);CHKERRQ(ierr); 136ddc54837SHong Zhang ierr = VecScatterCreate(red->xdup,is1,x,is2,&red->scatterout);CHKERRQ(ierr); 137fcfd50ebSBarry Smith ierr = ISDestroy(&is1);CHKERRQ(ierr); 138fcfd50ebSBarry Smith ierr = ISDestroy(&is2);CHKERRQ(ierr); 1391d79065fSBarry Smith ierr = PetscFree2(idx1,idx2);CHKERRQ(ierr); 140ddc54837SHong Zhang ierr = VecDestroy(&x);CHKERRQ(ierr); 1411b81debcSHong Zhang } 142ab661555SHong Zhang } else { /* !red->useparallelmat */ 143ab661555SHong Zhang ierr = KSPSetOperators(red->ksp,pc->mat,pc->pmat,pc->flag);CHKERRQ(ierr); 1441b81debcSHong Zhang } 145ab661555SHong Zhang } else { /* pc->setupcalled */ 1464b9ad928SBarry Smith if (red->useparallelmat) { 147ab661555SHong Zhang MatReuse reuse; 1481b81debcSHong Zhang /* grab the parallel matrix and put it into processors of a subcomminicator */ 1491b81debcSHong Zhang /*--------------------------------------------------------------------------*/ 150ab661555SHong Zhang if (pc->flag == DIFFERENT_NONZERO_PATTERN) { 1514b9ad928SBarry Smith /* destroy old matrices */ 1526bf464f9SBarry Smith ierr = MatDestroy(&red->pmats);CHKERRQ(ierr); 153ab661555SHong Zhang reuse = MAT_INITIAL_MATRIX; 1544b9ad928SBarry Smith } else { 155ab661555SHong Zhang reuse = MAT_REUSE_MATRIX; 156ab661555SHong Zhang } 157ab661555SHong Zhang ierr = MatGetRedundantMatrix_MPIAIJ_interlaced(pc->pmat,red->psubcomm->n,red->psubcomm,reuse,&red->pmats);CHKERRQ(ierr); 158ab661555SHong Zhang ierr = KSPSetOperators(red->ksp,red->pmats,red->pmats,pc->flag);CHKERRQ(ierr); 159ab661555SHong Zhang } else { /* !red->useparallelmat */ 16090f1c854SHong Zhang ierr = KSPSetOperators(red->ksp,pc->mat,pc->pmat,pc->flag);CHKERRQ(ierr); 1614b9ad928SBarry Smith } 162ab661555SHong Zhang } 1631b81debcSHong Zhang 1640c24e6a1SHong Zhang if (pc->setfromoptionscalled) { 1653e065800SHong Zhang ierr = KSPSetFromOptions(red->ksp);CHKERRQ(ierr); 1660c24e6a1SHong Zhang } 1673e065800SHong Zhang ierr = KSPSetUp(red->ksp);CHKERRQ(ierr); 1684b9ad928SBarry Smith PetscFunctionReturn(0); 1694b9ad928SBarry Smith } 1704b9ad928SBarry Smith 1714b9ad928SBarry Smith #undef __FUNCT__ 1724b9ad928SBarry Smith #define __FUNCT__ "PCApply_Redundant" 1736849ba73SBarry Smith static PetscErrorCode PCApply_Redundant(PC pc,Vec x,Vec y) 1744b9ad928SBarry Smith { 1754b9ad928SBarry Smith PC_Redundant *red = (PC_Redundant*)pc->data; 176dfbe8321SBarry Smith PetscErrorCode ierr; 1773f457be1SHong Zhang PetscScalar *array; 1784b9ad928SBarry Smith 1794b9ad928SBarry Smith PetscFunctionBegin; 180ddc54837SHong Zhang if (!red->useparallelmat) { 181ddc54837SHong Zhang ierr = KSPSolve(red->ksp,x,y);CHKERRQ(ierr); 182ddc54837SHong Zhang PetscFunctionReturn(0); 183ddc54837SHong Zhang } 184ddc54837SHong Zhang 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; 2181b81debcSHong Zhang if (red->useparallelmat) { 2196bf464f9SBarry Smith ierr = VecScatterDestroy(&red->scatterin);CHKERRQ(ierr); 2206bf464f9SBarry Smith ierr = VecScatterDestroy(&red->scatterout);CHKERRQ(ierr); 2216bf464f9SBarry Smith ierr = VecDestroy(&red->ysub);CHKERRQ(ierr); 2226bf464f9SBarry Smith ierr = VecDestroy(&red->xsub);CHKERRQ(ierr); 2236bf464f9SBarry Smith ierr = VecDestroy(&red->xdup);CHKERRQ(ierr); 2246bf464f9SBarry Smith ierr = VecDestroy(&red->ydup);CHKERRQ(ierr); 2251b81debcSHong Zhang } 2266bf464f9SBarry Smith ierr = MatDestroy(&red->pmats);CHKERRQ(ierr); 2271b81debcSHong Zhang ierr = KSPReset(red->ksp);CHKERRQ(ierr); 2281ea5a559SBarry Smith PetscFunctionReturn(0); 2291ea5a559SBarry Smith } 2301ea5a559SBarry Smith 2311ea5a559SBarry Smith #undef __FUNCT__ 2321ea5a559SBarry Smith #define __FUNCT__ "PCDestroy_Redundant" 2331ea5a559SBarry Smith static PetscErrorCode PCDestroy_Redundant(PC pc) 2341ea5a559SBarry Smith { 2351ea5a559SBarry Smith PC_Redundant *red = (PC_Redundant*)pc->data; 2361ea5a559SBarry Smith PetscErrorCode ierr; 2371ea5a559SBarry Smith 2381ea5a559SBarry Smith PetscFunctionBegin; 2391ea5a559SBarry Smith ierr = PCReset_Redundant(pc);CHKERRQ(ierr); 2406bf464f9SBarry Smith ierr = KSPDestroy(&red->ksp);CHKERRQ(ierr); 2416bf464f9SBarry Smith ierr = PetscSubcommDestroy(&red->psubcomm);CHKERRQ(ierr); 242c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 2434b9ad928SBarry Smith PetscFunctionReturn(0); 2444b9ad928SBarry Smith } 2454b9ad928SBarry Smith 2464b9ad928SBarry Smith #undef __FUNCT__ 2474b9ad928SBarry Smith #define __FUNCT__ "PCSetFromOptions_Redundant" 2486849ba73SBarry Smith static PetscErrorCode PCSetFromOptions_Redundant(PC pc) 2494b9ad928SBarry Smith { 250a98ce0f4SHong Zhang PetscErrorCode ierr; 251a98ce0f4SHong Zhang PC_Redundant *red = (PC_Redundant*)pc->data; 252a98ce0f4SHong Zhang 2534b9ad928SBarry Smith PetscFunctionBegin; 254a98ce0f4SHong Zhang ierr = PetscOptionsHead("Redundant options");CHKERRQ(ierr); 25509a6bc64SHong Zhang ierr = PetscOptionsInt("-pc_redundant_number","Number of redundant pc","PCRedundantSetNumber",red->nsubcomm,&red->nsubcomm,0);CHKERRQ(ierr); 256a98ce0f4SHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 2574b9ad928SBarry Smith PetscFunctionReturn(0); 2584b9ad928SBarry Smith } 2594b9ad928SBarry Smith 2604b9ad928SBarry Smith #undef __FUNCT__ 26109a6bc64SHong Zhang #define __FUNCT__ "PCRedundantSetNumber_Redundant" 262f7a08781SBarry Smith static PetscErrorCode PCRedundantSetNumber_Redundant(PC pc,PetscInt nreds) 26309a6bc64SHong Zhang { 26409a6bc64SHong Zhang PC_Redundant *red = (PC_Redundant*)pc->data; 26509a6bc64SHong Zhang 26609a6bc64SHong Zhang PetscFunctionBegin; 26709a6bc64SHong Zhang red->nsubcomm = nreds; 26809a6bc64SHong Zhang PetscFunctionReturn(0); 26909a6bc64SHong Zhang } 27009a6bc64SHong Zhang 27109a6bc64SHong Zhang #undef __FUNCT__ 27209a6bc64SHong Zhang #define __FUNCT__ "PCRedundantSetNumber" 27309a6bc64SHong Zhang /*@ 27409a6bc64SHong Zhang PCRedundantSetNumber - Sets the number of redundant preconditioner contexts. 27509a6bc64SHong Zhang 2763f9fe445SBarry Smith Logically Collective on PC 27709a6bc64SHong Zhang 27809a6bc64SHong Zhang Input Parameters: 27909a6bc64SHong Zhang + pc - the preconditioner context 2809b21d695SBarry Smith - nredundant - number of redundant preconditioner contexts; for example if you are using 64 MPI processes and 2819b21d695SBarry Smith use an nredundant of 4 there will be 4 parallel solves each on 16 = 64/4 processes. 28209a6bc64SHong Zhang 28309a6bc64SHong Zhang Level: advanced 28409a6bc64SHong Zhang 28509a6bc64SHong Zhang .keywords: PC, redundant solve 28609a6bc64SHong Zhang @*/ 2877087cfbeSBarry Smith PetscErrorCode PCRedundantSetNumber(PC pc,PetscInt nredundant) 28809a6bc64SHong Zhang { 2894ac538c5SBarry Smith PetscErrorCode ierr; 29009a6bc64SHong Zhang 29109a6bc64SHong Zhang PetscFunctionBegin; 2920700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 293ce94432eSBarry Smith if (nredundant <= 0) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG, "num of redundant pc %D must be positive",nredundant); 2944ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCRedundantSetNumber_C",(PC,PetscInt),(pc,nredundant));CHKERRQ(ierr); 29509a6bc64SHong Zhang PetscFunctionReturn(0); 29609a6bc64SHong Zhang } 29709a6bc64SHong Zhang 29809a6bc64SHong Zhang #undef __FUNCT__ 2994b9ad928SBarry Smith #define __FUNCT__ "PCRedundantSetScatter_Redundant" 300f7a08781SBarry Smith static PetscErrorCode PCRedundantSetScatter_Redundant(PC pc,VecScatter in,VecScatter out) 3014b9ad928SBarry Smith { 3024b9ad928SBarry Smith PC_Redundant *red = (PC_Redundant*)pc->data; 303dfbe8321SBarry Smith PetscErrorCode ierr; 3044b9ad928SBarry Smith 3054b9ad928SBarry Smith PetscFunctionBegin; 3064b9ad928SBarry Smith ierr = PetscObjectReference((PetscObject)in);CHKERRQ(ierr); 3076bf464f9SBarry Smith ierr = VecScatterDestroy(&red->scatterin);CHKERRQ(ierr); 3082fa5cd67SKarl Rupp 309c3122656SLisandro Dalcin red->scatterin = in; 3102fa5cd67SKarl Rupp 3114b9ad928SBarry Smith ierr = PetscObjectReference((PetscObject)out);CHKERRQ(ierr); 3126bf464f9SBarry Smith ierr = VecScatterDestroy(&red->scatterout);CHKERRQ(ierr); 313c3122656SLisandro Dalcin red->scatterout = out; 3144b9ad928SBarry Smith PetscFunctionReturn(0); 3154b9ad928SBarry Smith } 3164b9ad928SBarry Smith 3174b9ad928SBarry Smith #undef __FUNCT__ 3184b9ad928SBarry Smith #define __FUNCT__ "PCRedundantSetScatter" 3194b9ad928SBarry Smith /*@ 3204b9ad928SBarry Smith PCRedundantSetScatter - Sets the scatter used to copy values into the 3214b9ad928SBarry Smith redundant local solve and the scatter to move them back into the global 3224b9ad928SBarry Smith vector. 3234b9ad928SBarry Smith 3243f9fe445SBarry Smith Logically Collective on PC 3254b9ad928SBarry Smith 3264b9ad928SBarry Smith Input Parameters: 3274b9ad928SBarry Smith + pc - the preconditioner context 3284b9ad928SBarry Smith . in - the scatter to move the values in 3294b9ad928SBarry Smith - out - the scatter to move them out 3304b9ad928SBarry Smith 3314b9ad928SBarry Smith Level: advanced 3324b9ad928SBarry Smith 3334b9ad928SBarry Smith .keywords: PC, redundant solve 3344b9ad928SBarry Smith @*/ 3357087cfbeSBarry Smith PetscErrorCode PCRedundantSetScatter(PC pc,VecScatter in,VecScatter out) 3364b9ad928SBarry Smith { 3374ac538c5SBarry Smith PetscErrorCode ierr; 3384b9ad928SBarry Smith 3394b9ad928SBarry Smith PetscFunctionBegin; 3400700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 3410700a824SBarry Smith PetscValidHeaderSpecific(in,VEC_SCATTER_CLASSID,2); 3420700a824SBarry Smith PetscValidHeaderSpecific(out,VEC_SCATTER_CLASSID,3); 3434ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCRedundantSetScatter_C",(PC,VecScatter,VecScatter),(pc,in,out));CHKERRQ(ierr); 3444b9ad928SBarry Smith PetscFunctionReturn(0); 3454b9ad928SBarry Smith } 3464b9ad928SBarry Smith 3474b9ad928SBarry Smith #undef __FUNCT__ 34883ab6a24SBarry Smith #define __FUNCT__ "PCRedundantGetKSP_Redundant" 349f7a08781SBarry Smith static PetscErrorCode PCRedundantGetKSP_Redundant(PC pc,KSP *innerksp) 3504b9ad928SBarry Smith { 3515f06b7aaSBarry Smith PetscErrorCode ierr; 3524b9ad928SBarry Smith PC_Redundant *red = (PC_Redundant*)pc->data; 3535f06b7aaSBarry Smith MPI_Comm comm,subcomm; 3545f06b7aaSBarry Smith const char *prefix; 3554b9ad928SBarry Smith 3564b9ad928SBarry Smith PetscFunctionBegin; 3575f06b7aaSBarry Smith if (!red->psubcomm) { 3585f06b7aaSBarry Smith ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 359d8a68f86SHong Zhang ierr = PetscSubcommCreate(comm,&red->psubcomm);CHKERRQ(ierr); 360d8a68f86SHong Zhang ierr = PetscSubcommSetNumber(red->psubcomm,red->nsubcomm);CHKERRQ(ierr); 361d8a68f86SHong Zhang ierr = PetscSubcommSetType(red->psubcomm,PETSC_SUBCOMM_INTERLACED);CHKERRQ(ierr); 3625f06b7aaSBarry Smith ierr = PetscLogObjectMemory(pc,sizeof(PetscSubcomm));CHKERRQ(ierr); 3635f06b7aaSBarry Smith 3645f06b7aaSBarry Smith /* create a new PC that processors in each subcomm have copy of */ 3655f06b7aaSBarry Smith subcomm = red->psubcomm->comm; 3662fa5cd67SKarl Rupp 3675f06b7aaSBarry Smith ierr = KSPCreate(subcomm,&red->ksp);CHKERRQ(ierr); 3685f06b7aaSBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)red->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 3695f06b7aaSBarry Smith ierr = PetscLogObjectParent(pc,red->ksp);CHKERRQ(ierr); 3705f06b7aaSBarry Smith ierr = KSPSetType(red->ksp,KSPPREONLY);CHKERRQ(ierr); 3715f06b7aaSBarry Smith ierr = KSPGetPC(red->ksp,&red->pc);CHKERRQ(ierr); 3725f06b7aaSBarry Smith ierr = PCSetType(red->pc,PCLU);CHKERRQ(ierr); 3735f06b7aaSBarry Smith 3745f06b7aaSBarry Smith ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 3755f06b7aaSBarry Smith ierr = KSPSetOptionsPrefix(red->ksp,prefix);CHKERRQ(ierr); 3765f06b7aaSBarry Smith ierr = KSPAppendOptionsPrefix(red->ksp,"redundant_");CHKERRQ(ierr); 3775f06b7aaSBarry Smith } 37883ab6a24SBarry Smith *innerksp = red->ksp; 3794b9ad928SBarry Smith PetscFunctionReturn(0); 3804b9ad928SBarry Smith } 3814b9ad928SBarry Smith 3824b9ad928SBarry Smith #undef __FUNCT__ 38383ab6a24SBarry Smith #define __FUNCT__ "PCRedundantGetKSP" 3844b9ad928SBarry Smith /*@ 38583ab6a24SBarry Smith PCRedundantGetKSP - Gets the less parallel KSP 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: 39383ab6a24SBarry Smith . innerksp - the KSP on the smaller set of processes 3944b9ad928SBarry Smith 3954b9ad928SBarry Smith Level: advanced 3964b9ad928SBarry Smith 3974b9ad928SBarry Smith .keywords: PC, redundant solve 3984b9ad928SBarry Smith @*/ 39983ab6a24SBarry Smith PetscErrorCode PCRedundantGetKSP(PC pc,KSP *innerksp) 4004b9ad928SBarry Smith { 4014ac538c5SBarry Smith PetscErrorCode ierr; 4024b9ad928SBarry Smith 4034b9ad928SBarry Smith PetscFunctionBegin; 4040700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 40583ab6a24SBarry Smith PetscValidPointer(innerksp,2); 40683ab6a24SBarry Smith ierr = PetscTryMethod(pc,"PCRedundantGetKSP_C",(PC,KSP*),(pc,innerksp));CHKERRQ(ierr); 4074b9ad928SBarry Smith PetscFunctionReturn(0); 4084b9ad928SBarry Smith } 4094b9ad928SBarry Smith 4104b9ad928SBarry Smith #undef __FUNCT__ 4114b9ad928SBarry Smith #define __FUNCT__ "PCRedundantGetOperators_Redundant" 412f7a08781SBarry Smith static PetscErrorCode PCRedundantGetOperators_Redundant(PC pc,Mat *mat,Mat *pmat) 4134b9ad928SBarry Smith { 4144b9ad928SBarry Smith PC_Redundant *red = (PC_Redundant*)pc->data; 4154b9ad928SBarry Smith 4164b9ad928SBarry Smith PetscFunctionBegin; 417b3804887SHong Zhang if (mat) *mat = red->pmats; 418b3804887SHong Zhang if (pmat) *pmat = red->pmats; 4194b9ad928SBarry Smith PetscFunctionReturn(0); 4204b9ad928SBarry Smith } 4214b9ad928SBarry Smith 4224b9ad928SBarry Smith #undef __FUNCT__ 4234b9ad928SBarry Smith #define __FUNCT__ "PCRedundantGetOperators" 4244b9ad928SBarry Smith /*@ 4254b9ad928SBarry Smith PCRedundantGetOperators - gets the sequential matrix and preconditioner matrix 4264b9ad928SBarry Smith 4274b9ad928SBarry Smith Not Collective 4284b9ad928SBarry Smith 4294b9ad928SBarry Smith Input Parameter: 4304b9ad928SBarry Smith . pc - the preconditioner context 4314b9ad928SBarry Smith 4324b9ad928SBarry Smith Output Parameters: 4334b9ad928SBarry Smith + mat - the matrix 4344b9ad928SBarry Smith - pmat - the (possibly different) preconditioner matrix 4354b9ad928SBarry Smith 4364b9ad928SBarry Smith Level: advanced 4374b9ad928SBarry Smith 4384b9ad928SBarry Smith .keywords: PC, redundant solve 4394b9ad928SBarry Smith @*/ 4407087cfbeSBarry Smith PetscErrorCode PCRedundantGetOperators(PC pc,Mat *mat,Mat *pmat) 4414b9ad928SBarry Smith { 4424ac538c5SBarry Smith PetscErrorCode ierr; 4434b9ad928SBarry Smith 4444b9ad928SBarry Smith PetscFunctionBegin; 4450700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 4464482741eSBarry Smith if (mat) PetscValidPointer(mat,2); 4474482741eSBarry Smith if (pmat) PetscValidPointer(pmat,3); 4484ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCRedundantGetOperators_C",(PC,Mat*,Mat*),(pc,mat,pmat));CHKERRQ(ierr); 4494b9ad928SBarry Smith PetscFunctionReturn(0); 4504b9ad928SBarry Smith } 4514b9ad928SBarry Smith 4524b9ad928SBarry Smith /* -------------------------------------------------------------------------------------*/ 45337a17b4dSBarry Smith /*MC 45483ab6a24SBarry Smith PCREDUNDANT - Runs a KSP solver with preconditioner for the entire problem on subgroups of processors 45537a17b4dSBarry Smith 45683ab6a24SBarry Smith Options for the redundant preconditioners can be set with -redundant_pc_xxx for the redundant KSP with -redundant_ksp_xxx 45737a17b4dSBarry Smith 45809391456SBarry Smith Options Database: 4599b21d695SBarry Smith . -pc_redundant_number <n> - number of redundant solves, for example if you are using 64 MPI processes and 4609b21d695SBarry Smith use an n of 4 there will be 4 parallel solves each on 16 = 64/4 processes. 46109391456SBarry Smith 46237a17b4dSBarry Smith Level: intermediate 46337a17b4dSBarry Smith 46483ab6a24SBarry Smith Notes: The default KSP is preonly and the default PC is LU. 46583ab6a24SBarry Smith 46683ab6a24SBarry Smith Developer Notes: Note that PCSetInitialGuessNonzero() is not used by this class but likely should be. 4679cfaa89bSBarry Smith 46837a17b4dSBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PCRedundantSetScatter(), 46983ab6a24SBarry Smith PCRedundantGetKSP(), PCRedundantGetOperators(), PCRedundantSetNumber() 47037a17b4dSBarry Smith M*/ 47137a17b4dSBarry Smith 4724b9ad928SBarry Smith #undef __FUNCT__ 4734b9ad928SBarry Smith #define __FUNCT__ "PCCreate_Redundant" 4748cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_Redundant(PC pc) 4754b9ad928SBarry Smith { 476dfbe8321SBarry Smith PetscErrorCode ierr; 4774b9ad928SBarry Smith PC_Redundant *red; 47869db28dcSHong Zhang PetscMPIInt size; 4793f457be1SHong Zhang 4804b9ad928SBarry Smith PetscFunctionBegin; 48138f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_Redundant,&red);CHKERRQ(ierr); 482ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr); 4832fa5cd67SKarl Rupp 48469db28dcSHong Zhang red->nsubcomm = size; 4854b9ad928SBarry Smith red->useparallelmat = PETSC_TRUE; 4861fbd8f88SHong Zhang pc->data = (void*)red; 4874b9ad928SBarry Smith 4884b9ad928SBarry Smith pc->ops->apply = PCApply_Redundant; 4894b9ad928SBarry Smith pc->ops->applytranspose = 0; 4904b9ad928SBarry Smith pc->ops->setup = PCSetUp_Redundant; 4914b9ad928SBarry Smith pc->ops->destroy = PCDestroy_Redundant; 4921ea5a559SBarry Smith pc->ops->reset = PCReset_Redundant; 4934b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_Redundant; 4944b9ad928SBarry Smith pc->ops->view = PCView_Redundant; 4952fa5cd67SKarl Rupp 496bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCRedundantSetScatter_C",PCRedundantSetScatter_Redundant);CHKERRQ(ierr); 497bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCRedundantSetNumber_C",PCRedundantSetNumber_Redundant);CHKERRQ(ierr); 498bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCRedundantGetKSP_C",PCRedundantGetKSP_Redundant);CHKERRQ(ierr); 499bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCRedundantGetOperators_C",PCRedundantGetOperators_Redundant);CHKERRQ(ierr); 5004b9ad928SBarry Smith PetscFunctionReturn(0); 5014b9ad928SBarry Smith } 502b2573a8aSBarry Smith 503