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; 2713f74950SBarry Smith PetscMPIInt rank; 2832077d6dSBarry Smith PetscTruth iascii,isstring; 29a47c9f9aSHong Zhang PetscViewer sviewer,subviewer; 301fbd8f88SHong Zhang PetscInt color = red->psubcomm->color; 314b9ad928SBarry Smith 324b9ad928SBarry Smith PetscFunctionBegin; 337adad957SLisandro Dalcin ierr = MPI_Comm_rank(((PetscObject)pc)->comm,&rank);CHKERRQ(ierr); 3432077d6dSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 354b9ad928SBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr); 3632077d6dSBarry Smith if (iascii) { 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); 39a98ce0f4SHong Zhang if (!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); 45a98ce0f4SHong Zhang } else if (isstring) { /* not test it yet! */ 464b9ad928SBarry Smith ierr = PetscViewerStringSPrintf(viewer," Redundant solver preconditioner");CHKERRQ(ierr); 474b9ad928SBarry Smith ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr); 484b9ad928SBarry Smith if (!rank) { 493e065800SHong Zhang ierr = KSPView(red->ksp,sviewer);CHKERRQ(ierr); 504b9ad928SBarry Smith } 514b9ad928SBarry Smith ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr); 524b9ad928SBarry Smith } else { 5379a5c55eSBarry Smith SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PC redundant",((PetscObject)viewer)->type_name); 544b9ad928SBarry Smith } 554b9ad928SBarry Smith PetscFunctionReturn(0); 564b9ad928SBarry Smith } 574b9ad928SBarry Smith 587c4f633dSBarry Smith #include "private/matimpl.h" /*I "petscmat.h" I*/ 594b9ad928SBarry Smith #undef __FUNCT__ 604b9ad928SBarry Smith #define __FUNCT__ "PCSetUp_Redundant" 616849ba73SBarry Smith static PetscErrorCode PCSetUp_Redundant(PC pc) 624b9ad928SBarry Smith { 634b9ad928SBarry Smith PC_Redundant *red = (PC_Redundant*)pc->data; 64dfbe8321SBarry Smith PetscErrorCode ierr; 6545fc02eaSBarry Smith PetscInt mstart,mend,mlocal,m,mlocal_sub,rstart_sub,rend_sub,mloc_sub; 6613f74950SBarry Smith PetscMPIInt size; 674b9ad928SBarry Smith MatReuse reuse = MAT_INITIAL_MATRIX; 684b9ad928SBarry Smith MatStructure str = DIFFERENT_NONZERO_PATTERN; 697adad957SLisandro Dalcin MPI_Comm comm = ((PetscObject)pc)->comm,subcomm; 7023ce1328SBarry Smith Vec vec; 713f457be1SHong Zhang PetscMPIInt subsize,subrank; 721fbd8f88SHong Zhang const char *prefix; 733f457be1SHong Zhang 744b9ad928SBarry Smith PetscFunctionBegin; 7523ce1328SBarry Smith ierr = MatGetVecs(pc->pmat,&vec,0);CHKERRQ(ierr); 7623ce1328SBarry Smith ierr = VecGetSize(vec,&m);CHKERRQ(ierr); 771fbd8f88SHong Zhang 784b9ad928SBarry Smith if (!pc->setupcalled) { 79*5f06b7aaSBarry Smith if (!red->psubcomm) { 801fbd8f88SHong Zhang ierr = PetscSubcommCreate(comm,red->nsubcomm,&red->psubcomm);CHKERRQ(ierr); 811fbd8f88SHong Zhang ierr = PetscLogObjectMemory(pc,sizeof(PetscSubcomm));CHKERRQ(ierr); 821fbd8f88SHong Zhang 831fbd8f88SHong Zhang /* create a new PC that processors in each subcomm have copy of */ 840d7810c8SBarry Smith subcomm = red->psubcomm->comm; 85*5f06b7aaSBarry Smith ierr = KSPCreate(subcomm,&red->ksp);CHKERRQ(ierr); 86*5f06b7aaSBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)red->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 87*5f06b7aaSBarry Smith ierr = PetscLogObjectParent(pc,red->ksp);CHKERRQ(ierr); 88*5f06b7aaSBarry Smith ierr = KSPSetType(red->ksp,KSPPREONLY);CHKERRQ(ierr); 89*5f06b7aaSBarry Smith ierr = KSPGetPC(red->ksp,&red->pc);CHKERRQ(ierr); 90cf52b8b1SHong Zhang ierr = PCSetType(red->pc,PCLU);CHKERRQ(ierr); 91cf52b8b1SHong Zhang 921fbd8f88SHong Zhang ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 93*5f06b7aaSBarry Smith ierr = KSPSetOptionsPrefix(red->ksp,prefix);CHKERRQ(ierr); 94*5f06b7aaSBarry Smith ierr = KSPAppendOptionsPrefix(red->ksp,"redundant_");CHKERRQ(ierr); 95*5f06b7aaSBarry Smith } else { 96*5f06b7aaSBarry Smith subcomm = red->psubcomm->comm; 97*5f06b7aaSBarry Smith } 981fbd8f88SHong Zhang 993f457be1SHong Zhang /* create working vectors xsub/ysub and xdup/ydup */ 10023ce1328SBarry Smith ierr = VecGetLocalSize(vec,&mlocal);CHKERRQ(ierr); 1013f457be1SHong Zhang ierr = VecGetOwnershipRange(vec,&mstart,&mend);CHKERRQ(ierr); 1024b9ad928SBarry Smith 1033f457be1SHong Zhang /* get local size of xsub/ysub */ 1041fbd8f88SHong Zhang ierr = MPI_Comm_size(subcomm,&subsize);CHKERRQ(ierr); 1051fbd8f88SHong Zhang ierr = MPI_Comm_rank(subcomm,&subrank);CHKERRQ(ierr); 106d0f46423SBarry Smith rstart_sub = pc->pmat->rmap->range[red->psubcomm->n*subrank]; /* rstart in xsub/ysub */ 1073f457be1SHong Zhang if (subrank+1 < subsize){ 108d0f46423SBarry Smith rend_sub = pc->pmat->rmap->range[red->psubcomm->n*(subrank+1)]; 1093f457be1SHong Zhang } else { 1103f457be1SHong Zhang rend_sub = m; 1113f457be1SHong Zhang } 1123f457be1SHong Zhang mloc_sub = rend_sub - rstart_sub; 1131fbd8f88SHong Zhang ierr = VecCreateMPI(subcomm,mloc_sub,PETSC_DECIDE,&red->ysub);CHKERRQ(ierr); 1143f457be1SHong Zhang /* create xsub with empty local arrays, because xdup's arrays will be placed into it */ 1151fbd8f88SHong Zhang ierr = VecCreateMPIWithArray(subcomm,mloc_sub,PETSC_DECIDE,PETSC_NULL,&red->xsub);CHKERRQ(ierr); 1163f457be1SHong Zhang 1173f457be1SHong Zhang /* create xdup and ydup. ydup has empty local arrays because ysub's arrays will be place into it. 1187adad957SLisandro Dalcin Note: we use communicator dupcomm, not ((PetscObject)pc)->comm! */ 1191fbd8f88SHong Zhang ierr = VecCreateMPI(red->psubcomm->dupparent,mloc_sub,PETSC_DECIDE,&red->xdup);CHKERRQ(ierr); 1201fbd8f88SHong Zhang ierr = VecCreateMPIWithArray(red->psubcomm->dupparent,mloc_sub,PETSC_DECIDE,PETSC_NULL,&red->ydup);CHKERRQ(ierr); 1213f457be1SHong Zhang 1223f457be1SHong Zhang /* create vec scatters */ 1233f457be1SHong Zhang if (!red->scatterin){ 1243f457be1SHong Zhang IS is1,is2; 1253f457be1SHong Zhang PetscInt *idx1,*idx2,i,j,k; 12645fc02eaSBarry Smith 1271fbd8f88SHong Zhang ierr = PetscMalloc(2*red->psubcomm->n*mlocal*sizeof(PetscInt),&idx1);CHKERRQ(ierr); 1281fbd8f88SHong Zhang idx2 = idx1 + red->psubcomm->n*mlocal; 1293f457be1SHong Zhang j = 0; 1301fbd8f88SHong Zhang for (k=0; k<red->psubcomm->n; k++){ 1313f457be1SHong Zhang for (i=mstart; i<mend; i++){ 1323f457be1SHong Zhang idx1[j] = i; 1333f457be1SHong Zhang idx2[j++] = i + m*k; 1343f457be1SHong Zhang } 1353f457be1SHong Zhang } 1361fbd8f88SHong Zhang ierr = ISCreateGeneral(comm,red->psubcomm->n*mlocal,idx1,&is1);CHKERRQ(ierr); 1371fbd8f88SHong Zhang ierr = ISCreateGeneral(comm,red->psubcomm->n*mlocal,idx2,&is2);CHKERRQ(ierr); 1383f457be1SHong Zhang ierr = VecScatterCreate(vec,is1,red->xdup,is2,&red->scatterin);CHKERRQ(ierr); 1393f457be1SHong Zhang ierr = ISDestroy(is1);CHKERRQ(ierr); 1403f457be1SHong Zhang ierr = ISDestroy(is2);CHKERRQ(ierr); 1413f457be1SHong Zhang 1421fbd8f88SHong Zhang ierr = ISCreateStride(comm,mlocal,mstart+ red->psubcomm->color*m,1,&is1);CHKERRQ(ierr); 1433f457be1SHong Zhang ierr = ISCreateStride(comm,mlocal,mstart,1,&is2);CHKERRQ(ierr); 1443f457be1SHong Zhang ierr = VecScatterCreate(red->xdup,is1,vec,is2,&red->scatterout);CHKERRQ(ierr); 1453f457be1SHong Zhang ierr = ISDestroy(is1);CHKERRQ(ierr); 1463f457be1SHong Zhang ierr = ISDestroy(is2);CHKERRQ(ierr); 1473f457be1SHong Zhang ierr = PetscFree(idx1);CHKERRQ(ierr); 1484b9ad928SBarry Smith } 1494b9ad928SBarry Smith } 15023ce1328SBarry Smith ierr = VecDestroy(vec);CHKERRQ(ierr); 1514b9ad928SBarry Smith 1524b9ad928SBarry Smith /* if pmatrix set by user is sequential then we do not need to gather the parallel matrix */ 1533f457be1SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1544b9ad928SBarry Smith if (size == 1) { 1554b9ad928SBarry Smith red->useparallelmat = PETSC_FALSE; 1564b9ad928SBarry Smith } 1574b9ad928SBarry Smith 1584b9ad928SBarry Smith if (red->useparallelmat) { 1594b9ad928SBarry Smith if (pc->setupcalled == 1 && pc->flag == DIFFERENT_NONZERO_PATTERN) { 1604b9ad928SBarry Smith /* destroy old matrices */ 1614b9ad928SBarry Smith if (red->pmats) { 162b3804887SHong Zhang ierr = MatDestroy(red->pmats);CHKERRQ(ierr); 1634b9ad928SBarry Smith } 1644b9ad928SBarry Smith } else if (pc->setupcalled == 1) { 1654b9ad928SBarry Smith reuse = MAT_REUSE_MATRIX; 1664b9ad928SBarry Smith str = SAME_NONZERO_PATTERN; 1674b9ad928SBarry Smith } 1684b9ad928SBarry Smith 1693f457be1SHong Zhang /* grab the parallel matrix and put it into processors of a subcomminicator */ 170f664ae05SHong Zhang /*--------------------------------------------------------------------------*/ 171f664ae05SHong Zhang ierr = VecGetLocalSize(red->ysub,&mlocal_sub);CHKERRQ(ierr); 17269db28dcSHong Zhang ierr = MatGetRedundantMatrix(pc->pmat,red->psubcomm->n,red->psubcomm->comm,mlocal_sub,reuse,&red->pmats);CHKERRQ(ierr); 1733f457be1SHong Zhang /* tell PC of the subcommunicator its operators */ 17490f1c854SHong Zhang ierr = KSPSetOperators(red->ksp,red->pmats,red->pmats,str);CHKERRQ(ierr); 1754b9ad928SBarry Smith } else { 17690f1c854SHong Zhang ierr = KSPSetOperators(red->ksp,pc->mat,pc->pmat,pc->flag);CHKERRQ(ierr); 1774b9ad928SBarry Smith } 1780c24e6a1SHong Zhang if (pc->setfromoptionscalled){ 1793e065800SHong Zhang ierr = KSPSetFromOptions(red->ksp);CHKERRQ(ierr); 1800c24e6a1SHong Zhang } 1813e065800SHong Zhang ierr = KSPSetUp(red->ksp);CHKERRQ(ierr); 1824b9ad928SBarry Smith PetscFunctionReturn(0); 1834b9ad928SBarry Smith } 1844b9ad928SBarry Smith 1854b9ad928SBarry Smith #undef __FUNCT__ 1864b9ad928SBarry Smith #define __FUNCT__ "PCApply_Redundant" 1876849ba73SBarry Smith static PetscErrorCode PCApply_Redundant(PC pc,Vec x,Vec y) 1884b9ad928SBarry Smith { 1894b9ad928SBarry Smith PC_Redundant *red = (PC_Redundant*)pc->data; 190dfbe8321SBarry Smith PetscErrorCode ierr; 1913f457be1SHong Zhang PetscScalar *array; 1924b9ad928SBarry Smith 1934b9ad928SBarry Smith PetscFunctionBegin; 1943f457be1SHong Zhang /* scatter x to xdup */ 195ca9f406cSSatish Balay ierr = VecScatterBegin(red->scatterin,x,red->xdup,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 196ca9f406cSSatish Balay ierr = VecScatterEnd(red->scatterin,x,red->xdup,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1973f457be1SHong Zhang 1983f457be1SHong Zhang /* place xdup's local array into xsub */ 1993f457be1SHong Zhang ierr = VecGetArray(red->xdup,&array);CHKERRQ(ierr); 2003f457be1SHong Zhang ierr = VecPlaceArray(red->xsub,(const PetscScalar*)array);CHKERRQ(ierr); 2014b9ad928SBarry Smith 2024b9ad928SBarry Smith /* apply preconditioner on each processor */ 2033f457be1SHong Zhang ierr = PCApply(red->pc,red->xsub,red->ysub);CHKERRQ(ierr); 2043f457be1SHong Zhang ierr = VecResetArray(red->xsub);CHKERRQ(ierr); 2053f457be1SHong Zhang ierr = VecRestoreArray(red->xdup,&array);CHKERRQ(ierr); 2064b9ad928SBarry Smith 2073f457be1SHong Zhang /* place ysub's local array into ydup */ 2083f457be1SHong Zhang ierr = VecGetArray(red->ysub,&array);CHKERRQ(ierr); 2093f457be1SHong Zhang ierr = VecPlaceArray(red->ydup,(const PetscScalar*)array);CHKERRQ(ierr); 2103f457be1SHong Zhang 2113f457be1SHong Zhang /* scatter ydup to y */ 212ca9f406cSSatish Balay ierr = VecScatterBegin(red->scatterout,red->ydup,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 213ca9f406cSSatish Balay ierr = VecScatterEnd(red->scatterout,red->ydup,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 2143f457be1SHong Zhang ierr = VecResetArray(red->ydup);CHKERRQ(ierr); 2153f457be1SHong Zhang ierr = VecRestoreArray(red->ysub,&array);CHKERRQ(ierr); 2164b9ad928SBarry Smith PetscFunctionReturn(0); 2174b9ad928SBarry Smith } 2184b9ad928SBarry Smith 2194b9ad928SBarry Smith #undef __FUNCT__ 2204b9ad928SBarry Smith #define __FUNCT__ "PCDestroy_Redundant" 2216849ba73SBarry Smith static PetscErrorCode PCDestroy_Redundant(PC pc) 2224b9ad928SBarry Smith { 2234b9ad928SBarry Smith PC_Redundant *red = (PC_Redundant*)pc->data; 224dfbe8321SBarry Smith PetscErrorCode ierr; 2254b9ad928SBarry Smith 2264b9ad928SBarry Smith PetscFunctionBegin; 2274b9ad928SBarry Smith if (red->scatterin) {ierr = VecScatterDestroy(red->scatterin);CHKERRQ(ierr);} 2284b9ad928SBarry Smith if (red->scatterout) {ierr = VecScatterDestroy(red->scatterout);CHKERRQ(ierr);} 2293f457be1SHong Zhang if (red->ysub) {ierr = VecDestroy(red->ysub);CHKERRQ(ierr);} 2303f457be1SHong Zhang if (red->xsub) {ierr = VecDestroy(red->xsub);CHKERRQ(ierr);} 2313f457be1SHong Zhang if (red->xdup) {ierr = VecDestroy(red->xdup);CHKERRQ(ierr);} 2323f457be1SHong Zhang if (red->ydup) {ierr = VecDestroy(red->ydup);CHKERRQ(ierr);} 233b3804887SHong Zhang if (red->pmats) { 234b3804887SHong Zhang ierr = MatDestroy(red->pmats);CHKERRQ(ierr); 2353f457be1SHong Zhang } 236e5a9bf91SBarry Smith if (red->psubcomm) {ierr = PetscSubcommDestroy(red->psubcomm);CHKERRQ(ierr);} 2373e065800SHong Zhang if (red->ksp) {ierr = KSPDestroy(red->ksp);CHKERRQ(ierr);} 2384b9ad928SBarry Smith ierr = PetscFree(red);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 EXTERN_C_BEGIN 2574b9ad928SBarry Smith #undef __FUNCT__ 25809a6bc64SHong Zhang #define __FUNCT__ "PCRedundantSetNumber_Redundant" 25909a6bc64SHong Zhang PetscErrorCode PETSCKSP_DLLEXPORT PCRedundantSetNumber_Redundant(PC pc,PetscInt nreds) 26009a6bc64SHong Zhang { 26109a6bc64SHong Zhang PC_Redundant *red = (PC_Redundant*)pc->data; 26209a6bc64SHong Zhang 26309a6bc64SHong Zhang PetscFunctionBegin; 26409a6bc64SHong Zhang red->nsubcomm = nreds; 26509a6bc64SHong Zhang PetscFunctionReturn(0); 26609a6bc64SHong Zhang } 26709a6bc64SHong Zhang EXTERN_C_END 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 27409a6bc64SHong Zhang 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 @*/ 28509a6bc64SHong Zhang PetscErrorCode PETSCKSP_DLLEXPORT PCRedundantSetNumber(PC pc,PetscInt nredundant) 28609a6bc64SHong Zhang { 28709a6bc64SHong Zhang PetscErrorCode ierr,(*f)(PC,PetscInt); 28809a6bc64SHong Zhang 28909a6bc64SHong Zhang PetscFunctionBegin; 29009a6bc64SHong Zhang PetscValidHeaderSpecific(pc,PC_COOKIE,1); 29109a6bc64SHong Zhang if (nredundant <= 0) SETERRQ1(PETSC_ERR_ARG_WRONG, "num of redundant pc %D must be positive",nredundant); 29209a6bc64SHong Zhang ierr = PetscObjectQueryFunction((PetscObject)pc,"PCRedundantSetNumber_C",(void (**)(void))&f);CHKERRQ(ierr); 29309a6bc64SHong Zhang if (f) { 29409a6bc64SHong Zhang ierr = (*f)(pc,nredundant);CHKERRQ(ierr); 29509a6bc64SHong Zhang } 29609a6bc64SHong Zhang PetscFunctionReturn(0); 29709a6bc64SHong Zhang } 29809a6bc64SHong Zhang 29909a6bc64SHong Zhang EXTERN_C_BEGIN 30009a6bc64SHong Zhang #undef __FUNCT__ 3014b9ad928SBarry Smith #define __FUNCT__ "PCRedundantSetScatter_Redundant" 302dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCRedundantSetScatter_Redundant(PC pc,VecScatter in,VecScatter out) 3034b9ad928SBarry Smith { 3044b9ad928SBarry Smith PC_Redundant *red = (PC_Redundant*)pc->data; 305dfbe8321SBarry Smith PetscErrorCode ierr; 3064b9ad928SBarry Smith 3074b9ad928SBarry Smith PetscFunctionBegin; 3084b9ad928SBarry Smith ierr = PetscObjectReference((PetscObject)in);CHKERRQ(ierr); 309c3122656SLisandro Dalcin if (red->scatterin) { ierr = VecScatterDestroy(red->scatterin);CHKERRQ(ierr); } 310c3122656SLisandro Dalcin red->scatterin = in; 3114b9ad928SBarry Smith ierr = PetscObjectReference((PetscObject)out);CHKERRQ(ierr); 312c3122656SLisandro Dalcin if (red->scatterout) { ierr = VecScatterDestroy(red->scatterout);CHKERRQ(ierr); } 313c3122656SLisandro Dalcin red->scatterout = out; 3144b9ad928SBarry Smith PetscFunctionReturn(0); 3154b9ad928SBarry Smith } 3164b9ad928SBarry Smith EXTERN_C_END 3174b9ad928SBarry Smith 3184b9ad928SBarry Smith #undef __FUNCT__ 3194b9ad928SBarry Smith #define __FUNCT__ "PCRedundantSetScatter" 3204b9ad928SBarry Smith /*@ 3214b9ad928SBarry Smith PCRedundantSetScatter - Sets the scatter used to copy values into the 3224b9ad928SBarry Smith redundant local solve and the scatter to move them back into the global 3234b9ad928SBarry Smith vector. 3244b9ad928SBarry Smith 3254b9ad928SBarry Smith Collective on PC 3264b9ad928SBarry Smith 3274b9ad928SBarry Smith Input Parameters: 3284b9ad928SBarry Smith + pc - the preconditioner context 3294b9ad928SBarry Smith . in - the scatter to move the values in 3304b9ad928SBarry Smith - out - the scatter to move them out 3314b9ad928SBarry Smith 3324b9ad928SBarry Smith Level: advanced 3334b9ad928SBarry Smith 3344b9ad928SBarry Smith .keywords: PC, redundant solve 3354b9ad928SBarry Smith @*/ 336dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCRedundantSetScatter(PC pc,VecScatter in,VecScatter out) 3374b9ad928SBarry Smith { 338dfbe8321SBarry Smith PetscErrorCode ierr,(*f)(PC,VecScatter,VecScatter); 3394b9ad928SBarry Smith 3404b9ad928SBarry Smith PetscFunctionBegin; 3414482741eSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 3424482741eSBarry Smith PetscValidHeaderSpecific(in,VEC_SCATTER_COOKIE,2); 3434482741eSBarry Smith PetscValidHeaderSpecific(out,VEC_SCATTER_COOKIE,3); 3444b9ad928SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCRedundantSetScatter_C",(void (**)(void))&f);CHKERRQ(ierr); 3454b9ad928SBarry Smith if (f) { 3464b9ad928SBarry Smith ierr = (*f)(pc,in,out);CHKERRQ(ierr); 3474b9ad928SBarry Smith } 3484b9ad928SBarry Smith PetscFunctionReturn(0); 3494b9ad928SBarry Smith } 3504b9ad928SBarry Smith 3514b9ad928SBarry Smith EXTERN_C_BEGIN 3524b9ad928SBarry Smith #undef __FUNCT__ 3534b9ad928SBarry Smith #define __FUNCT__ "PCRedundantGetPC_Redundant" 354dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCRedundantGetPC_Redundant(PC pc,PC *innerpc) 3554b9ad928SBarry Smith { 356*5f06b7aaSBarry Smith PetscErrorCode ierr; 3574b9ad928SBarry Smith PC_Redundant *red = (PC_Redundant*)pc->data; 358*5f06b7aaSBarry Smith MPI_Comm comm,subcomm; 359*5f06b7aaSBarry Smith const char *prefix; 3604b9ad928SBarry Smith 3614b9ad928SBarry Smith PetscFunctionBegin; 362*5f06b7aaSBarry Smith if (!red->psubcomm) { 363*5f06b7aaSBarry Smith ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 364*5f06b7aaSBarry Smith ierr = PetscSubcommCreate(comm,red->nsubcomm,&red->psubcomm);CHKERRQ(ierr); 365*5f06b7aaSBarry Smith ierr = PetscLogObjectMemory(pc,sizeof(PetscSubcomm));CHKERRQ(ierr); 366*5f06b7aaSBarry Smith 367*5f06b7aaSBarry Smith /* create a new PC that processors in each subcomm have copy of */ 368*5f06b7aaSBarry Smith subcomm = red->psubcomm->comm; 369*5f06b7aaSBarry Smith ierr = KSPCreate(subcomm,&red->ksp);CHKERRQ(ierr); 370*5f06b7aaSBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)red->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 371*5f06b7aaSBarry Smith ierr = PetscLogObjectParent(pc,red->ksp);CHKERRQ(ierr); 372*5f06b7aaSBarry Smith ierr = KSPSetType(red->ksp,KSPPREONLY);CHKERRQ(ierr); 373*5f06b7aaSBarry Smith ierr = KSPGetPC(red->ksp,&red->pc);CHKERRQ(ierr); 374*5f06b7aaSBarry Smith ierr = PCSetType(red->pc,PCLU);CHKERRQ(ierr); 375*5f06b7aaSBarry Smith 376*5f06b7aaSBarry Smith ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 377*5f06b7aaSBarry Smith ierr = KSPSetOptionsPrefix(red->ksp,prefix);CHKERRQ(ierr); 378*5f06b7aaSBarry Smith ierr = KSPAppendOptionsPrefix(red->ksp,"redundant_");CHKERRQ(ierr); 379*5f06b7aaSBarry Smith } 380*5f06b7aaSBarry Smith 381*5f06b7aaSBarry Smith ierr = KSPGetPC(red->ksp,innerpc);CHKERRQ(ierr); 3824b9ad928SBarry Smith PetscFunctionReturn(0); 3834b9ad928SBarry Smith } 3844b9ad928SBarry Smith EXTERN_C_END 3854b9ad928SBarry Smith 3864b9ad928SBarry Smith #undef __FUNCT__ 3874b9ad928SBarry Smith #define __FUNCT__ "PCRedundantGetPC" 3884b9ad928SBarry Smith /*@ 3894b9ad928SBarry Smith PCRedundantGetPC - Gets the sequential PC created by the redundant PC. 3904b9ad928SBarry Smith 3914b9ad928SBarry Smith Not Collective 3924b9ad928SBarry Smith 3934b9ad928SBarry Smith Input Parameter: 3944b9ad928SBarry Smith . pc - the preconditioner context 3954b9ad928SBarry Smith 3964b9ad928SBarry Smith Output Parameter: 3974b9ad928SBarry Smith . innerpc - the sequential PC 3984b9ad928SBarry Smith 3994b9ad928SBarry Smith Level: advanced 4004b9ad928SBarry Smith 4014b9ad928SBarry Smith .keywords: PC, redundant solve 4024b9ad928SBarry Smith @*/ 403dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCRedundantGetPC(PC pc,PC *innerpc) 4044b9ad928SBarry Smith { 405dfbe8321SBarry Smith PetscErrorCode ierr,(*f)(PC,PC*); 4064b9ad928SBarry Smith 4074b9ad928SBarry Smith PetscFunctionBegin; 4084482741eSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 4094482741eSBarry Smith PetscValidPointer(innerpc,2); 4104b9ad928SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCRedundantGetPC_C",(void (**)(void))&f);CHKERRQ(ierr); 4114b9ad928SBarry Smith if (f) { 4124b9ad928SBarry Smith ierr = (*f)(pc,innerpc);CHKERRQ(ierr); 4134b9ad928SBarry Smith } 4144b9ad928SBarry Smith PetscFunctionReturn(0); 4154b9ad928SBarry Smith } 4164b9ad928SBarry Smith 4174b9ad928SBarry Smith EXTERN_C_BEGIN 4184b9ad928SBarry Smith #undef __FUNCT__ 4194b9ad928SBarry Smith #define __FUNCT__ "PCRedundantGetOperators_Redundant" 420dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCRedundantGetOperators_Redundant(PC pc,Mat *mat,Mat *pmat) 4214b9ad928SBarry Smith { 4224b9ad928SBarry Smith PC_Redundant *red = (PC_Redundant*)pc->data; 4234b9ad928SBarry Smith 4244b9ad928SBarry Smith PetscFunctionBegin; 425b3804887SHong Zhang if (mat) *mat = red->pmats; 426b3804887SHong Zhang if (pmat) *pmat = red->pmats; 4274b9ad928SBarry Smith PetscFunctionReturn(0); 4284b9ad928SBarry Smith } 4294b9ad928SBarry Smith EXTERN_C_END 4304b9ad928SBarry Smith 4314b9ad928SBarry Smith #undef __FUNCT__ 4324b9ad928SBarry Smith #define __FUNCT__ "PCRedundantGetOperators" 4334b9ad928SBarry Smith /*@ 4344b9ad928SBarry Smith PCRedundantGetOperators - gets the sequential matrix and preconditioner matrix 4354b9ad928SBarry Smith 4364b9ad928SBarry Smith Not Collective 4374b9ad928SBarry Smith 4384b9ad928SBarry Smith Input Parameter: 4394b9ad928SBarry Smith . pc - the preconditioner context 4404b9ad928SBarry Smith 4414b9ad928SBarry Smith Output Parameters: 4424b9ad928SBarry Smith + mat - the matrix 4434b9ad928SBarry Smith - pmat - the (possibly different) preconditioner matrix 4444b9ad928SBarry Smith 4454b9ad928SBarry Smith Level: advanced 4464b9ad928SBarry Smith 4474b9ad928SBarry Smith .keywords: PC, redundant solve 4484b9ad928SBarry Smith @*/ 449dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCRedundantGetOperators(PC pc,Mat *mat,Mat *pmat) 4504b9ad928SBarry Smith { 451dfbe8321SBarry Smith PetscErrorCode ierr,(*f)(PC,Mat*,Mat*); 4524b9ad928SBarry Smith 4534b9ad928SBarry Smith PetscFunctionBegin; 4544482741eSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 4554482741eSBarry Smith if (mat) PetscValidPointer(mat,2); 4564482741eSBarry Smith if (pmat) PetscValidPointer(pmat,3); 4574b9ad928SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCRedundantGetOperators_C",(void (**)(void))&f);CHKERRQ(ierr); 4584b9ad928SBarry Smith if (f) { 4594b9ad928SBarry Smith ierr = (*f)(pc,mat,pmat);CHKERRQ(ierr); 4604b9ad928SBarry Smith } 4614b9ad928SBarry Smith PetscFunctionReturn(0); 4624b9ad928SBarry Smith } 4634b9ad928SBarry Smith 4644b9ad928SBarry Smith /* -------------------------------------------------------------------------------------*/ 46537a17b4dSBarry Smith /*MC 466a98ce0f4SHong Zhang PCREDUNDANT - Runs a preconditioner for the entire problem on subgroups of processors 46737a17b4dSBarry Smith 46837a17b4dSBarry Smith Options for the redundant preconditioners can be set with -redundant_pc_xxx 46937a17b4dSBarry Smith 47009391456SBarry Smith Options Database: 4719b21d695SBarry Smith . -pc_redundant_number <n> - number of redundant solves, for example if you are using 64 MPI processes and 4729b21d695SBarry Smith use an n of 4 there will be 4 parallel solves each on 16 = 64/4 processes. 47309391456SBarry Smith 47437a17b4dSBarry Smith Level: intermediate 47537a17b4dSBarry Smith 47637a17b4dSBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PCRedundantSetScatter(), 4779b21d695SBarry Smith PCRedundantGetPC(), PCRedundantGetOperators(), PCRedundantSetNumber() 47837a17b4dSBarry Smith M*/ 47937a17b4dSBarry Smith 4804b9ad928SBarry Smith EXTERN_C_BEGIN 4814b9ad928SBarry Smith #undef __FUNCT__ 4824b9ad928SBarry Smith #define __FUNCT__ "PCCreate_Redundant" 483dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_Redundant(PC pc) 4844b9ad928SBarry Smith { 485dfbe8321SBarry Smith PetscErrorCode ierr; 4864b9ad928SBarry Smith PC_Redundant *red; 48769db28dcSHong Zhang PetscMPIInt size; 4883f457be1SHong Zhang 4894b9ad928SBarry Smith PetscFunctionBegin; 49038f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_Redundant,&red);CHKERRQ(ierr); 4917adad957SLisandro Dalcin ierr = MPI_Comm_size(((PetscObject)pc)->comm,&size);CHKERRQ(ierr); 49269db28dcSHong Zhang red->nsubcomm = size; 4934b9ad928SBarry Smith red->useparallelmat = PETSC_TRUE; 4941fbd8f88SHong Zhang pc->data = (void*)red; 4954b9ad928SBarry Smith 4964b9ad928SBarry Smith pc->ops->apply = PCApply_Redundant; 4974b9ad928SBarry Smith pc->ops->applytranspose = 0; 4984b9ad928SBarry Smith pc->ops->setup = PCSetUp_Redundant; 4994b9ad928SBarry Smith pc->ops->destroy = PCDestroy_Redundant; 5004b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_Redundant; 5014b9ad928SBarry Smith pc->ops->view = PCView_Redundant; 5024b9ad928SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCRedundantSetScatter_C","PCRedundantSetScatter_Redundant", 5034b9ad928SBarry Smith PCRedundantSetScatter_Redundant);CHKERRQ(ierr); 50409a6bc64SHong Zhang ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCRedundantSetNumber_C","PCRedundantSetNumber_Redundant", 50509a6bc64SHong Zhang PCRedundantSetNumber_Redundant);CHKERRQ(ierr); 5064b9ad928SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCRedundantGetPC_C","PCRedundantGetPC_Redundant", 5074b9ad928SBarry Smith PCRedundantGetPC_Redundant);CHKERRQ(ierr); 5084b9ad928SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCRedundantGetOperators_C","PCRedundantGetOperators_Redundant", 5094b9ad928SBarry Smith PCRedundantGetOperators_Redundant);CHKERRQ(ierr); 5104b9ad928SBarry Smith PetscFunctionReturn(0); 5114b9ad928SBarry Smith } 5124b9ad928SBarry Smith EXTERN_C_END 513