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; 571b81debcSHong Zhang PetscInt mstart,mend,mlocal,M; 5813f74950SBarry Smith PetscMPIInt size; 59ce94432eSBarry Smith MPI_Comm comm,subcomm; 60ddc54837SHong Zhang Vec x; 611fbd8f88SHong Zhang const char *prefix; 623f457be1SHong Zhang 634b9ad928SBarry Smith PetscFunctionBegin; 64ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 65ddc54837SHong Zhang 66ddc54837SHong Zhang /* if pmatrix set by user is sequential then we do not need to gather the parallel matrix */ 67ddc54837SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 68ddc54837SHong Zhang if (size == 1) red->useparallelmat = PETSC_FALSE; 691fbd8f88SHong Zhang 704b9ad928SBarry Smith if (!pc->setupcalled) { 711b81debcSHong Zhang PetscInt mloc_sub; 725f06b7aaSBarry Smith if (!red->psubcomm) { 73d8a68f86SHong Zhang ierr = PetscSubcommCreate(comm,&red->psubcomm);CHKERRQ(ierr); 74d8a68f86SHong Zhang ierr = PetscSubcommSetNumber(red->psubcomm,red->nsubcomm);CHKERRQ(ierr); 75a23b1e67SHong Zhang ierr = PetscSubcommSetType(red->psubcomm,PETSC_SUBCOMM_CONTIGUOUS);CHKERRQ(ierr); 76a23b1e67SHong Zhang /* enable runtime switch of psubcomm type, e.g., '-psubcomm_type interlaced */ 77a23b1e67SHong Zhang ierr = PetscSubcommSetFromOptions(red->psubcomm);CHKERRQ(ierr); 783bb1ff40SBarry Smith ierr = PetscLogObjectMemory((PetscObject)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); 853bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)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 */ 99b2bf6370SHong Zhang ierr = MatGetRedundantMatrix(pc->pmat,red->psubcomm->n,subcomm,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 1058b96b0d2SHong Zhang /* create working vectors xdup and ydup. 1068b96b0d2SHong Zhang xdup concatenates all xsub's contigously to form a mpi vector over dupcomm (see PetscSubcommCreate_interlaced()) 1078b96b0d2SHong Zhang ydup concatenates all ysub and has empty local arrays because ysub's arrays will be place into it. 108ce94432eSBarry Smith Note: we use communicator dupcomm, not PetscObjectComm((PetscObject)pc)! */ 1091b81debcSHong Zhang ierr = MatGetLocalSize(red->pmats,&mloc_sub,NULL);CHKERRQ(ierr); 1101fbd8f88SHong Zhang ierr = VecCreateMPI(red->psubcomm->dupparent,mloc_sub,PETSC_DECIDE,&red->xdup);CHKERRQ(ierr); 1110298fd71SBarry Smith ierr = VecCreateMPIWithArray(red->psubcomm->dupparent,1,mloc_sub,PETSC_DECIDE,NULL,&red->ydup);CHKERRQ(ierr); 1123f457be1SHong Zhang 113f68be91cSHong Zhang /* create vecscatters */ 114f68be91cSHong Zhang if (!red->scatterin) { /* efficiency of scatterin is independent from psubcomm_type! */ 1153f457be1SHong Zhang IS is1,is2; 1163f457be1SHong Zhang PetscInt *idx1,*idx2,i,j,k; 11745fc02eaSBarry Smith 1181b81debcSHong Zhang ierr = MatGetVecs(pc->pmat,&x,0);CHKERRQ(ierr); 1191b81debcSHong Zhang ierr = VecGetSize(x,&M);CHKERRQ(ierr); 1201b81debcSHong Zhang ierr = VecGetOwnershipRange(x,&mstart,&mend);CHKERRQ(ierr); 1211b81debcSHong Zhang mlocal = mend - mstart; 122dcca6d9dSJed Brown ierr = PetscMalloc2(red->psubcomm->n*mlocal,&idx1,red->psubcomm->n*mlocal,&idx2);CHKERRQ(ierr); 1233f457be1SHong Zhang j = 0; 1241fbd8f88SHong Zhang for (k=0; k<red->psubcomm->n; k++) { 1253f457be1SHong Zhang for (i=mstart; i<mend; i++) { 1263f457be1SHong Zhang idx1[j] = i; 127ddc54837SHong Zhang idx2[j++] = i + M*k; 1283f457be1SHong Zhang } 1293f457be1SHong Zhang } 13070b3c8c7SBarry Smith ierr = ISCreateGeneral(comm,red->psubcomm->n*mlocal,idx1,PETSC_COPY_VALUES,&is1);CHKERRQ(ierr); 13170b3c8c7SBarry Smith ierr = ISCreateGeneral(comm,red->psubcomm->n*mlocal,idx2,PETSC_COPY_VALUES,&is2);CHKERRQ(ierr); 132ddc54837SHong Zhang ierr = VecScatterCreate(x,is1,red->xdup,is2,&red->scatterin);CHKERRQ(ierr); 133fcfd50ebSBarry Smith ierr = ISDestroy(&is1);CHKERRQ(ierr); 134fcfd50ebSBarry Smith ierr = ISDestroy(&is2);CHKERRQ(ierr); 1353f457be1SHong Zhang 1366909748bSHong Zhang /* Impl below is good for PETSC_SUBCOMM_INTERLACED (no inter-process communication) and PETSC_SUBCOMM_CONTIGUOUS (communication within subcomm) */ 137ddc54837SHong Zhang ierr = ISCreateStride(comm,mlocal,mstart+ red->psubcomm->color*M,1,&is1);CHKERRQ(ierr); 1383f457be1SHong Zhang ierr = ISCreateStride(comm,mlocal,mstart,1,&is2);CHKERRQ(ierr); 139ddc54837SHong Zhang ierr = VecScatterCreate(red->xdup,is1,x,is2,&red->scatterout);CHKERRQ(ierr); 140fcfd50ebSBarry Smith ierr = ISDestroy(&is1);CHKERRQ(ierr); 141fcfd50ebSBarry Smith ierr = ISDestroy(&is2);CHKERRQ(ierr); 1421d79065fSBarry Smith ierr = PetscFree2(idx1,idx2);CHKERRQ(ierr); 143ddc54837SHong Zhang ierr = VecDestroy(&x);CHKERRQ(ierr); 1441b81debcSHong Zhang } 145ab661555SHong Zhang } else { /* !red->useparallelmat */ 146ab661555SHong Zhang ierr = KSPSetOperators(red->ksp,pc->mat,pc->pmat,pc->flag);CHKERRQ(ierr); 1471b81debcSHong Zhang } 148ab661555SHong Zhang } else { /* pc->setupcalled */ 1494b9ad928SBarry Smith if (red->useparallelmat) { 150ab661555SHong Zhang MatReuse reuse; 1511b81debcSHong Zhang /* grab the parallel matrix and put it into processors of a subcomminicator */ 1521b81debcSHong Zhang /*--------------------------------------------------------------------------*/ 153ab661555SHong Zhang if (pc->flag == DIFFERENT_NONZERO_PATTERN) { 1544b9ad928SBarry Smith /* destroy old matrices */ 1556bf464f9SBarry Smith ierr = MatDestroy(&red->pmats);CHKERRQ(ierr); 156ab661555SHong Zhang reuse = MAT_INITIAL_MATRIX; 1574b9ad928SBarry Smith } else { 158ab661555SHong Zhang reuse = MAT_REUSE_MATRIX; 159ab661555SHong Zhang } 160b2bf6370SHong Zhang ierr = MatGetRedundantMatrix(pc->pmat,red->psubcomm->n,red->psubcomm->comm,reuse,&red->pmats);CHKERRQ(ierr); 161ab661555SHong Zhang ierr = KSPSetOperators(red->ksp,red->pmats,red->pmats,pc->flag);CHKERRQ(ierr); 162ab661555SHong Zhang } else { /* !red->useparallelmat */ 16390f1c854SHong Zhang ierr = KSPSetOperators(red->ksp,pc->mat,pc->pmat,pc->flag);CHKERRQ(ierr); 1644b9ad928SBarry Smith } 165ab661555SHong Zhang } 1661b81debcSHong Zhang 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; 183ddc54837SHong Zhang if (!red->useparallelmat) { 184ddc54837SHong Zhang ierr = KSPSolve(red->ksp,x,y);CHKERRQ(ierr); 185ddc54837SHong Zhang PetscFunctionReturn(0); 186ddc54837SHong Zhang } 187ddc54837SHong Zhang 1883f457be1SHong Zhang /* scatter x to xdup */ 189ca9f406cSSatish Balay ierr = VecScatterBegin(red->scatterin,x,red->xdup,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 190ca9f406cSSatish Balay ierr = VecScatterEnd(red->scatterin,x,red->xdup,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 1913f457be1SHong Zhang 1923f457be1SHong Zhang /* place xdup's local array into xsub */ 1933f457be1SHong Zhang ierr = VecGetArray(red->xdup,&array);CHKERRQ(ierr); 1943f457be1SHong Zhang ierr = VecPlaceArray(red->xsub,(const PetscScalar*)array);CHKERRQ(ierr); 1954b9ad928SBarry Smith 1964b9ad928SBarry Smith /* apply preconditioner on each processor */ 19783ab6a24SBarry Smith ierr = KSPSolve(red->ksp,red->xsub,red->ysub);CHKERRQ(ierr); 1983f457be1SHong Zhang ierr = VecResetArray(red->xsub);CHKERRQ(ierr); 1993f457be1SHong Zhang ierr = VecRestoreArray(red->xdup,&array);CHKERRQ(ierr); 2004b9ad928SBarry Smith 2013f457be1SHong Zhang /* place ysub's local array into ydup */ 2023f457be1SHong Zhang ierr = VecGetArray(red->ysub,&array);CHKERRQ(ierr); 2033f457be1SHong Zhang ierr = VecPlaceArray(red->ydup,(const PetscScalar*)array);CHKERRQ(ierr); 2043f457be1SHong Zhang 2053f457be1SHong Zhang /* scatter ydup to y */ 206ca9f406cSSatish Balay ierr = VecScatterBegin(red->scatterout,red->ydup,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 207ca9f406cSSatish Balay ierr = VecScatterEnd(red->scatterout,red->ydup,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 2083f457be1SHong Zhang ierr = VecResetArray(red->ydup);CHKERRQ(ierr); 2093f457be1SHong Zhang ierr = VecRestoreArray(red->ysub,&array);CHKERRQ(ierr); 2104b9ad928SBarry Smith PetscFunctionReturn(0); 2114b9ad928SBarry Smith } 2124b9ad928SBarry Smith 2134b9ad928SBarry Smith #undef __FUNCT__ 214*d88bfacbSStefano Zampini #define __FUNCT__ "PCApplyTranspose_Redundant" 215*d88bfacbSStefano Zampini static PetscErrorCode PCApplyTranspose_Redundant(PC pc,Vec x,Vec y) 216*d88bfacbSStefano Zampini { 217*d88bfacbSStefano Zampini PC_Redundant *red = (PC_Redundant*)pc->data; 218*d88bfacbSStefano Zampini PetscErrorCode ierr; 219*d88bfacbSStefano Zampini PetscScalar *array; 220*d88bfacbSStefano Zampini 221*d88bfacbSStefano Zampini PetscFunctionBegin; 222*d88bfacbSStefano Zampini if (!red->useparallelmat) { 223*d88bfacbSStefano Zampini ierr = KSPSolveTranspose(red->ksp,x,y);CHKERRQ(ierr); 224*d88bfacbSStefano Zampini PetscFunctionReturn(0); 225*d88bfacbSStefano Zampini } 226*d88bfacbSStefano Zampini 227*d88bfacbSStefano Zampini /* scatter x to xdup */ 228*d88bfacbSStefano Zampini ierr = VecScatterBegin(red->scatterin,x,red->xdup,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 229*d88bfacbSStefano Zampini ierr = VecScatterEnd(red->scatterin,x,red->xdup,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 230*d88bfacbSStefano Zampini 231*d88bfacbSStefano Zampini /* place xdup's local array into xsub */ 232*d88bfacbSStefano Zampini ierr = VecGetArray(red->xdup,&array);CHKERRQ(ierr); 233*d88bfacbSStefano Zampini ierr = VecPlaceArray(red->xsub,(const PetscScalar*)array);CHKERRQ(ierr); 234*d88bfacbSStefano Zampini 235*d88bfacbSStefano Zampini /* apply preconditioner on each processor */ 236*d88bfacbSStefano Zampini ierr = KSPSolveTranspose(red->ksp,red->xsub,red->ysub);CHKERRQ(ierr); 237*d88bfacbSStefano Zampini ierr = VecResetArray(red->xsub);CHKERRQ(ierr); 238*d88bfacbSStefano Zampini ierr = VecRestoreArray(red->xdup,&array);CHKERRQ(ierr); 239*d88bfacbSStefano Zampini 240*d88bfacbSStefano Zampini /* place ysub's local array into ydup */ 241*d88bfacbSStefano Zampini ierr = VecGetArray(red->ysub,&array);CHKERRQ(ierr); 242*d88bfacbSStefano Zampini ierr = VecPlaceArray(red->ydup,(const PetscScalar*)array);CHKERRQ(ierr); 243*d88bfacbSStefano Zampini 244*d88bfacbSStefano Zampini /* scatter ydup to y */ 245*d88bfacbSStefano Zampini ierr = VecScatterBegin(red->scatterout,red->ydup,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 246*d88bfacbSStefano Zampini ierr = VecScatterEnd(red->scatterout,red->ydup,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 247*d88bfacbSStefano Zampini ierr = VecResetArray(red->ydup);CHKERRQ(ierr); 248*d88bfacbSStefano Zampini ierr = VecRestoreArray(red->ysub,&array);CHKERRQ(ierr); 249*d88bfacbSStefano Zampini PetscFunctionReturn(0); 250*d88bfacbSStefano Zampini } 251*d88bfacbSStefano Zampini 252*d88bfacbSStefano Zampini #undef __FUNCT__ 2531ea5a559SBarry Smith #define __FUNCT__ "PCReset_Redundant" 2541ea5a559SBarry Smith static PetscErrorCode PCReset_Redundant(PC pc) 2554b9ad928SBarry Smith { 2564b9ad928SBarry Smith PC_Redundant *red = (PC_Redundant*)pc->data; 257dfbe8321SBarry Smith PetscErrorCode ierr; 2584b9ad928SBarry Smith 2594b9ad928SBarry Smith PetscFunctionBegin; 2601b81debcSHong Zhang if (red->useparallelmat) { 2616bf464f9SBarry Smith ierr = VecScatterDestroy(&red->scatterin);CHKERRQ(ierr); 2626bf464f9SBarry Smith ierr = VecScatterDestroy(&red->scatterout);CHKERRQ(ierr); 2636bf464f9SBarry Smith ierr = VecDestroy(&red->ysub);CHKERRQ(ierr); 2646bf464f9SBarry Smith ierr = VecDestroy(&red->xsub);CHKERRQ(ierr); 2656bf464f9SBarry Smith ierr = VecDestroy(&red->xdup);CHKERRQ(ierr); 2666bf464f9SBarry Smith ierr = VecDestroy(&red->ydup);CHKERRQ(ierr); 2671b81debcSHong Zhang } 2686bf464f9SBarry Smith ierr = MatDestroy(&red->pmats);CHKERRQ(ierr); 2691b81debcSHong Zhang ierr = KSPReset(red->ksp);CHKERRQ(ierr); 2701ea5a559SBarry Smith PetscFunctionReturn(0); 2711ea5a559SBarry Smith } 2721ea5a559SBarry Smith 2731ea5a559SBarry Smith #undef __FUNCT__ 2741ea5a559SBarry Smith #define __FUNCT__ "PCDestroy_Redundant" 2751ea5a559SBarry Smith static PetscErrorCode PCDestroy_Redundant(PC pc) 2761ea5a559SBarry Smith { 2771ea5a559SBarry Smith PC_Redundant *red = (PC_Redundant*)pc->data; 2781ea5a559SBarry Smith PetscErrorCode ierr; 2791ea5a559SBarry Smith 2801ea5a559SBarry Smith PetscFunctionBegin; 2811ea5a559SBarry Smith ierr = PCReset_Redundant(pc);CHKERRQ(ierr); 2826bf464f9SBarry Smith ierr = KSPDestroy(&red->ksp);CHKERRQ(ierr); 2836bf464f9SBarry Smith ierr = PetscSubcommDestroy(&red->psubcomm);CHKERRQ(ierr); 284c31cb41cSBarry Smith ierr = PetscFree(pc->data);CHKERRQ(ierr); 2854b9ad928SBarry Smith PetscFunctionReturn(0); 2864b9ad928SBarry Smith } 2874b9ad928SBarry Smith 2884b9ad928SBarry Smith #undef __FUNCT__ 2894b9ad928SBarry Smith #define __FUNCT__ "PCSetFromOptions_Redundant" 2906849ba73SBarry Smith static PetscErrorCode PCSetFromOptions_Redundant(PC pc) 2914b9ad928SBarry Smith { 292a98ce0f4SHong Zhang PetscErrorCode ierr; 293a98ce0f4SHong Zhang PC_Redundant *red = (PC_Redundant*)pc->data; 294a98ce0f4SHong Zhang 2954b9ad928SBarry Smith PetscFunctionBegin; 296a98ce0f4SHong Zhang ierr = PetscOptionsHead("Redundant options");CHKERRQ(ierr); 29709a6bc64SHong Zhang ierr = PetscOptionsInt("-pc_redundant_number","Number of redundant pc","PCRedundantSetNumber",red->nsubcomm,&red->nsubcomm,0);CHKERRQ(ierr); 298a98ce0f4SHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 2994b9ad928SBarry Smith PetscFunctionReturn(0); 3004b9ad928SBarry Smith } 3014b9ad928SBarry Smith 3024b9ad928SBarry Smith #undef __FUNCT__ 30309a6bc64SHong Zhang #define __FUNCT__ "PCRedundantSetNumber_Redundant" 304f7a08781SBarry Smith static PetscErrorCode PCRedundantSetNumber_Redundant(PC pc,PetscInt nreds) 30509a6bc64SHong Zhang { 30609a6bc64SHong Zhang PC_Redundant *red = (PC_Redundant*)pc->data; 30709a6bc64SHong Zhang 30809a6bc64SHong Zhang PetscFunctionBegin; 30909a6bc64SHong Zhang red->nsubcomm = nreds; 31009a6bc64SHong Zhang PetscFunctionReturn(0); 31109a6bc64SHong Zhang } 31209a6bc64SHong Zhang 31309a6bc64SHong Zhang #undef __FUNCT__ 31409a6bc64SHong Zhang #define __FUNCT__ "PCRedundantSetNumber" 31509a6bc64SHong Zhang /*@ 31609a6bc64SHong Zhang PCRedundantSetNumber - Sets the number of redundant preconditioner contexts. 31709a6bc64SHong Zhang 3183f9fe445SBarry Smith Logically Collective on PC 31909a6bc64SHong Zhang 32009a6bc64SHong Zhang Input Parameters: 32109a6bc64SHong Zhang + pc - the preconditioner context 3229b21d695SBarry Smith - nredundant - number of redundant preconditioner contexts; for example if you are using 64 MPI processes and 3239b21d695SBarry Smith use an nredundant of 4 there will be 4 parallel solves each on 16 = 64/4 processes. 32409a6bc64SHong Zhang 32509a6bc64SHong Zhang Level: advanced 32609a6bc64SHong Zhang 32709a6bc64SHong Zhang .keywords: PC, redundant solve 32809a6bc64SHong Zhang @*/ 3297087cfbeSBarry Smith PetscErrorCode PCRedundantSetNumber(PC pc,PetscInt nredundant) 33009a6bc64SHong Zhang { 3314ac538c5SBarry Smith PetscErrorCode ierr; 33209a6bc64SHong Zhang 33309a6bc64SHong Zhang PetscFunctionBegin; 3340700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 335ce94432eSBarry Smith if (nredundant <= 0) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG, "num of redundant pc %D must be positive",nredundant); 3364ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCRedundantSetNumber_C",(PC,PetscInt),(pc,nredundant));CHKERRQ(ierr); 33709a6bc64SHong Zhang PetscFunctionReturn(0); 33809a6bc64SHong Zhang } 33909a6bc64SHong Zhang 34009a6bc64SHong Zhang #undef __FUNCT__ 3414b9ad928SBarry Smith #define __FUNCT__ "PCRedundantSetScatter_Redundant" 342f7a08781SBarry Smith static PetscErrorCode PCRedundantSetScatter_Redundant(PC pc,VecScatter in,VecScatter out) 3434b9ad928SBarry Smith { 3444b9ad928SBarry Smith PC_Redundant *red = (PC_Redundant*)pc->data; 345dfbe8321SBarry Smith PetscErrorCode ierr; 3464b9ad928SBarry Smith 3474b9ad928SBarry Smith PetscFunctionBegin; 3484b9ad928SBarry Smith ierr = PetscObjectReference((PetscObject)in);CHKERRQ(ierr); 3496bf464f9SBarry Smith ierr = VecScatterDestroy(&red->scatterin);CHKERRQ(ierr); 3502fa5cd67SKarl Rupp 351c3122656SLisandro Dalcin red->scatterin = in; 3522fa5cd67SKarl Rupp 3534b9ad928SBarry Smith ierr = PetscObjectReference((PetscObject)out);CHKERRQ(ierr); 3546bf464f9SBarry Smith ierr = VecScatterDestroy(&red->scatterout);CHKERRQ(ierr); 355c3122656SLisandro Dalcin red->scatterout = out; 3564b9ad928SBarry Smith PetscFunctionReturn(0); 3574b9ad928SBarry Smith } 3584b9ad928SBarry Smith 3594b9ad928SBarry Smith #undef __FUNCT__ 3604b9ad928SBarry Smith #define __FUNCT__ "PCRedundantSetScatter" 3614b9ad928SBarry Smith /*@ 3624b9ad928SBarry Smith PCRedundantSetScatter - Sets the scatter used to copy values into the 3634b9ad928SBarry Smith redundant local solve and the scatter to move them back into the global 3644b9ad928SBarry Smith vector. 3654b9ad928SBarry Smith 3663f9fe445SBarry Smith Logically Collective on PC 3674b9ad928SBarry Smith 3684b9ad928SBarry Smith Input Parameters: 3694b9ad928SBarry Smith + pc - the preconditioner context 3704b9ad928SBarry Smith . in - the scatter to move the values in 3714b9ad928SBarry Smith - out - the scatter to move them out 3724b9ad928SBarry Smith 3734b9ad928SBarry Smith Level: advanced 3744b9ad928SBarry Smith 3754b9ad928SBarry Smith .keywords: PC, redundant solve 3764b9ad928SBarry Smith @*/ 3777087cfbeSBarry Smith PetscErrorCode PCRedundantSetScatter(PC pc,VecScatter in,VecScatter out) 3784b9ad928SBarry Smith { 3794ac538c5SBarry Smith PetscErrorCode ierr; 3804b9ad928SBarry Smith 3814b9ad928SBarry Smith PetscFunctionBegin; 3820700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 3830700a824SBarry Smith PetscValidHeaderSpecific(in,VEC_SCATTER_CLASSID,2); 3840700a824SBarry Smith PetscValidHeaderSpecific(out,VEC_SCATTER_CLASSID,3); 3854ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCRedundantSetScatter_C",(PC,VecScatter,VecScatter),(pc,in,out));CHKERRQ(ierr); 3864b9ad928SBarry Smith PetscFunctionReturn(0); 3874b9ad928SBarry Smith } 3884b9ad928SBarry Smith 3894b9ad928SBarry Smith #undef __FUNCT__ 39083ab6a24SBarry Smith #define __FUNCT__ "PCRedundantGetKSP_Redundant" 391f7a08781SBarry Smith static PetscErrorCode PCRedundantGetKSP_Redundant(PC pc,KSP *innerksp) 3924b9ad928SBarry Smith { 3935f06b7aaSBarry Smith PetscErrorCode ierr; 3944b9ad928SBarry Smith PC_Redundant *red = (PC_Redundant*)pc->data; 3955f06b7aaSBarry Smith MPI_Comm comm,subcomm; 3965f06b7aaSBarry Smith const char *prefix; 3974b9ad928SBarry Smith 3984b9ad928SBarry Smith PetscFunctionBegin; 3995f06b7aaSBarry Smith if (!red->psubcomm) { 4005f06b7aaSBarry Smith ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr); 401d8a68f86SHong Zhang ierr = PetscSubcommCreate(comm,&red->psubcomm);CHKERRQ(ierr); 402d8a68f86SHong Zhang ierr = PetscSubcommSetNumber(red->psubcomm,red->nsubcomm);CHKERRQ(ierr); 403d8a68f86SHong Zhang ierr = PetscSubcommSetType(red->psubcomm,PETSC_SUBCOMM_INTERLACED);CHKERRQ(ierr); 4043bb1ff40SBarry Smith ierr = PetscLogObjectMemory((PetscObject)pc,sizeof(PetscSubcomm));CHKERRQ(ierr); 4055f06b7aaSBarry Smith 4065f06b7aaSBarry Smith /* create a new PC that processors in each subcomm have copy of */ 4075f06b7aaSBarry Smith subcomm = red->psubcomm->comm; 4082fa5cd67SKarl Rupp 4095f06b7aaSBarry Smith ierr = KSPCreate(subcomm,&red->ksp);CHKERRQ(ierr); 4105f06b7aaSBarry Smith ierr = PetscObjectIncrementTabLevel((PetscObject)red->ksp,(PetscObject)pc,1);CHKERRQ(ierr); 4113bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)red->ksp);CHKERRQ(ierr); 4125f06b7aaSBarry Smith ierr = KSPSetType(red->ksp,KSPPREONLY);CHKERRQ(ierr); 4135f06b7aaSBarry Smith ierr = KSPGetPC(red->ksp,&red->pc);CHKERRQ(ierr); 4145f06b7aaSBarry Smith ierr = PCSetType(red->pc,PCLU);CHKERRQ(ierr); 4155f06b7aaSBarry Smith 4165f06b7aaSBarry Smith ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 4175f06b7aaSBarry Smith ierr = KSPSetOptionsPrefix(red->ksp,prefix);CHKERRQ(ierr); 4185f06b7aaSBarry Smith ierr = KSPAppendOptionsPrefix(red->ksp,"redundant_");CHKERRQ(ierr); 4195f06b7aaSBarry Smith } 42083ab6a24SBarry Smith *innerksp = red->ksp; 4214b9ad928SBarry Smith PetscFunctionReturn(0); 4224b9ad928SBarry Smith } 4234b9ad928SBarry Smith 4244b9ad928SBarry Smith #undef __FUNCT__ 42583ab6a24SBarry Smith #define __FUNCT__ "PCRedundantGetKSP" 4264b9ad928SBarry Smith /*@ 42783ab6a24SBarry Smith PCRedundantGetKSP - Gets the less parallel KSP created by the redundant PC. 4284b9ad928SBarry Smith 4294b9ad928SBarry Smith Not Collective 4304b9ad928SBarry Smith 4314b9ad928SBarry Smith Input Parameter: 4324b9ad928SBarry Smith . pc - the preconditioner context 4334b9ad928SBarry Smith 4344b9ad928SBarry Smith Output Parameter: 43583ab6a24SBarry Smith . innerksp - the KSP on the smaller set of processes 4364b9ad928SBarry Smith 4374b9ad928SBarry Smith Level: advanced 4384b9ad928SBarry Smith 4394b9ad928SBarry Smith .keywords: PC, redundant solve 4404b9ad928SBarry Smith @*/ 44183ab6a24SBarry Smith PetscErrorCode PCRedundantGetKSP(PC pc,KSP *innerksp) 4424b9ad928SBarry Smith { 4434ac538c5SBarry Smith PetscErrorCode ierr; 4444b9ad928SBarry Smith 4454b9ad928SBarry Smith PetscFunctionBegin; 4460700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 44783ab6a24SBarry Smith PetscValidPointer(innerksp,2); 44883ab6a24SBarry Smith ierr = PetscTryMethod(pc,"PCRedundantGetKSP_C",(PC,KSP*),(pc,innerksp));CHKERRQ(ierr); 4494b9ad928SBarry Smith PetscFunctionReturn(0); 4504b9ad928SBarry Smith } 4514b9ad928SBarry Smith 4524b9ad928SBarry Smith #undef __FUNCT__ 4534b9ad928SBarry Smith #define __FUNCT__ "PCRedundantGetOperators_Redundant" 454f7a08781SBarry Smith static PetscErrorCode PCRedundantGetOperators_Redundant(PC pc,Mat *mat,Mat *pmat) 4554b9ad928SBarry Smith { 4564b9ad928SBarry Smith PC_Redundant *red = (PC_Redundant*)pc->data; 4574b9ad928SBarry Smith 4584b9ad928SBarry Smith PetscFunctionBegin; 459b3804887SHong Zhang if (mat) *mat = red->pmats; 460b3804887SHong Zhang if (pmat) *pmat = red->pmats; 4614b9ad928SBarry Smith PetscFunctionReturn(0); 4624b9ad928SBarry Smith } 4634b9ad928SBarry Smith 4644b9ad928SBarry Smith #undef __FUNCT__ 4654b9ad928SBarry Smith #define __FUNCT__ "PCRedundantGetOperators" 4664b9ad928SBarry Smith /*@ 4674b9ad928SBarry Smith PCRedundantGetOperators - gets the sequential matrix and preconditioner matrix 4684b9ad928SBarry Smith 4694b9ad928SBarry Smith Not Collective 4704b9ad928SBarry Smith 4714b9ad928SBarry Smith Input Parameter: 4724b9ad928SBarry Smith . pc - the preconditioner context 4734b9ad928SBarry Smith 4744b9ad928SBarry Smith Output Parameters: 4754b9ad928SBarry Smith + mat - the matrix 4764b9ad928SBarry Smith - pmat - the (possibly different) preconditioner matrix 4774b9ad928SBarry Smith 4784b9ad928SBarry Smith Level: advanced 4794b9ad928SBarry Smith 4804b9ad928SBarry Smith .keywords: PC, redundant solve 4814b9ad928SBarry Smith @*/ 4827087cfbeSBarry Smith PetscErrorCode PCRedundantGetOperators(PC pc,Mat *mat,Mat *pmat) 4834b9ad928SBarry Smith { 4844ac538c5SBarry Smith PetscErrorCode ierr; 4854b9ad928SBarry Smith 4864b9ad928SBarry Smith PetscFunctionBegin; 4870700a824SBarry Smith PetscValidHeaderSpecific(pc,PC_CLASSID,1); 4884482741eSBarry Smith if (mat) PetscValidPointer(mat,2); 4894482741eSBarry Smith if (pmat) PetscValidPointer(pmat,3); 4904ac538c5SBarry Smith ierr = PetscTryMethod(pc,"PCRedundantGetOperators_C",(PC,Mat*,Mat*),(pc,mat,pmat));CHKERRQ(ierr); 4914b9ad928SBarry Smith PetscFunctionReturn(0); 4924b9ad928SBarry Smith } 4934b9ad928SBarry Smith 4944b9ad928SBarry Smith /* -------------------------------------------------------------------------------------*/ 49537a17b4dSBarry Smith /*MC 49683ab6a24SBarry Smith PCREDUNDANT - Runs a KSP solver with preconditioner for the entire problem on subgroups of processors 49737a17b4dSBarry Smith 49883ab6a24SBarry Smith Options for the redundant preconditioners can be set with -redundant_pc_xxx for the redundant KSP with -redundant_ksp_xxx 49937a17b4dSBarry Smith 50009391456SBarry Smith Options Database: 5019b21d695SBarry Smith . -pc_redundant_number <n> - number of redundant solves, for example if you are using 64 MPI processes and 5029b21d695SBarry Smith use an n of 4 there will be 4 parallel solves each on 16 = 64/4 processes. 50309391456SBarry Smith 50437a17b4dSBarry Smith Level: intermediate 50537a17b4dSBarry Smith 50683ab6a24SBarry Smith Notes: The default KSP is preonly and the default PC is LU. 50783ab6a24SBarry Smith 50883ab6a24SBarry Smith Developer Notes: Note that PCSetInitialGuessNonzero() is not used by this class but likely should be. 5099cfaa89bSBarry Smith 51037a17b4dSBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PCRedundantSetScatter(), 51183ab6a24SBarry Smith PCRedundantGetKSP(), PCRedundantGetOperators(), PCRedundantSetNumber() 51237a17b4dSBarry Smith M*/ 51337a17b4dSBarry Smith 5144b9ad928SBarry Smith #undef __FUNCT__ 5154b9ad928SBarry Smith #define __FUNCT__ "PCCreate_Redundant" 5168cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_Redundant(PC pc) 5174b9ad928SBarry Smith { 518dfbe8321SBarry Smith PetscErrorCode ierr; 5194b9ad928SBarry Smith PC_Redundant *red; 52069db28dcSHong Zhang PetscMPIInt size; 5213f457be1SHong Zhang 5224b9ad928SBarry Smith PetscFunctionBegin; 523b00a9115SJed Brown ierr = PetscNewLog(pc,&red);CHKERRQ(ierr); 524ce94432eSBarry Smith ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr); 5252fa5cd67SKarl Rupp 52669db28dcSHong Zhang red->nsubcomm = size; 5274b9ad928SBarry Smith red->useparallelmat = PETSC_TRUE; 5281fbd8f88SHong Zhang pc->data = (void*)red; 5294b9ad928SBarry Smith 5304b9ad928SBarry Smith pc->ops->apply = PCApply_Redundant; 531*d88bfacbSStefano Zampini pc->ops->applytranspose = PCApplyTranspose_Redundant; 5324b9ad928SBarry Smith pc->ops->setup = PCSetUp_Redundant; 5334b9ad928SBarry Smith pc->ops->destroy = PCDestroy_Redundant; 5341ea5a559SBarry Smith pc->ops->reset = PCReset_Redundant; 5354b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_Redundant; 5364b9ad928SBarry Smith pc->ops->view = PCView_Redundant; 5372fa5cd67SKarl Rupp 538bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCRedundantSetScatter_C",PCRedundantSetScatter_Redundant);CHKERRQ(ierr); 539bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCRedundantSetNumber_C",PCRedundantSetNumber_Redundant);CHKERRQ(ierr); 540bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCRedundantGetKSP_C",PCRedundantGetKSP_Redundant);CHKERRQ(ierr); 541bdf89e91SBarry Smith ierr = PetscObjectComposeFunction((PetscObject)pc,"PCRedundantGetOperators_C",PCRedundantGetOperators_Redundant);CHKERRQ(ierr); 5424b9ad928SBarry Smith PetscFunctionReturn(0); 5434b9ad928SBarry Smith } 544b2573a8aSBarry Smith 545