xref: /petsc/src/ksp/pc/impls/redundant/redundant.c (revision 2fa5cd679192b9b390e47ae2d0650965e6b1d9fa)
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 */
117adad957SLisandro Dalcin   Vec          xsub,ysub;            /* vectors of a subcommunicator to hold parallel vectors of ((PetscObject)pc)->comm */
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;
617adad957SLisandro Dalcin   MPI_Comm       comm  = ((PetscObject)pc)->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;
6823ce1328SBarry Smith   ierr = MatGetVecs(pc->pmat,&vec,0);CHKERRQ(ierr);
6923ce1328SBarry Smith   ierr = VecGetSize(vec,&m);CHKERRQ(ierr);
701fbd8f88SHong Zhang 
714b9ad928SBarry Smith   if (!pc->setupcalled) {
725f06b7aaSBarry Smith     if (!red->psubcomm) {
73d8a68f86SHong Zhang       ierr = PetscSubcommCreate(comm,&red->psubcomm);CHKERRQ(ierr);
74d8a68f86SHong Zhang       ierr = PetscSubcommSetNumber(red->psubcomm,red->nsubcomm);CHKERRQ(ierr);
75d8a68f86SHong Zhang       ierr = PetscSubcommSetType(red->psubcomm,PETSC_SUBCOMM_INTERLACED);CHKERRQ(ierr);
761fbd8f88SHong Zhang       ierr = PetscLogObjectMemory(pc,sizeof(PetscSubcomm));CHKERRQ(ierr);
771fbd8f88SHong Zhang 
781fbd8f88SHong Zhang       /* create a new PC that processors in each subcomm have copy of */
790d7810c8SBarry Smith       subcomm = red->psubcomm->comm;
80*2fa5cd67SKarl Rupp 
815f06b7aaSBarry Smith       ierr = KSPCreate(subcomm,&red->ksp);CHKERRQ(ierr);
825f06b7aaSBarry Smith       ierr = PetscObjectIncrementTabLevel((PetscObject)red->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
835f06b7aaSBarry Smith       ierr = PetscLogObjectParent(pc,red->ksp);CHKERRQ(ierr);
845f06b7aaSBarry Smith       ierr = KSPSetType(red->ksp,KSPPREONLY);CHKERRQ(ierr);
855f06b7aaSBarry Smith       ierr = KSPGetPC(red->ksp,&red->pc);CHKERRQ(ierr);
86cf52b8b1SHong Zhang       ierr = PCSetType(red->pc,PCLU);CHKERRQ(ierr);
87cf52b8b1SHong Zhang 
881fbd8f88SHong Zhang       ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
895f06b7aaSBarry Smith       ierr = KSPSetOptionsPrefix(red->ksp,prefix);CHKERRQ(ierr);
905f06b7aaSBarry Smith       ierr = KSPAppendOptionsPrefix(red->ksp,"redundant_");CHKERRQ(ierr);
91*2fa5cd67SKarl Rupp     } else subcomm = red->psubcomm->comm;
921fbd8f88SHong Zhang 
933f457be1SHong Zhang     /* create working vectors xsub/ysub and xdup/ydup */
9423ce1328SBarry Smith     ierr = VecGetLocalSize(vec,&mlocal);CHKERRQ(ierr);
953f457be1SHong Zhang     ierr = VecGetOwnershipRange(vec,&mstart,&mend);CHKERRQ(ierr);
964b9ad928SBarry Smith 
973f457be1SHong Zhang     /* get local size of xsub/ysub */
981fbd8f88SHong Zhang     ierr = MPI_Comm_size(subcomm,&subsize);CHKERRQ(ierr);
991fbd8f88SHong Zhang     ierr = MPI_Comm_rank(subcomm,&subrank);CHKERRQ(ierr);
100b862ddfaSBarry Smith     ierr = MatGetOwnershipRanges(pc->pmat,&range);CHKERRQ(ierr);
101b862ddfaSBarry Smith     rstart_sub = range[red->psubcomm->n*subrank]; /* rstart in xsub/ysub */
102*2fa5cd67SKarl Rupp     if (subrank+1 < subsize) rend_sub = range[red->psubcomm->n*(subrank+1)];
103*2fa5cd67SKarl Rupp     else rend_sub = m;
104*2fa5cd67SKarl Rupp 
1053f457be1SHong Zhang     mloc_sub = rend_sub - rstart_sub;
1061fbd8f88SHong Zhang     ierr     = VecCreateMPI(subcomm,mloc_sub,PETSC_DECIDE,&red->ysub);CHKERRQ(ierr);
1073f457be1SHong Zhang     /* create xsub with empty local arrays, because xdup's arrays will be placed into it */
108778a2246SBarry Smith     ierr = VecCreateMPIWithArray(subcomm,1,mloc_sub,PETSC_DECIDE,PETSC_NULL,&red->xsub);CHKERRQ(ierr);
1093f457be1SHong Zhang 
1103f457be1SHong Zhang     /* create xdup and ydup. ydup has empty local arrays because ysub's arrays will be place into it.
1117adad957SLisandro Dalcin        Note: we use communicator dupcomm, not ((PetscObject)pc)->comm! */
1121fbd8f88SHong Zhang     ierr = VecCreateMPI(red->psubcomm->dupparent,mloc_sub,PETSC_DECIDE,&red->xdup);CHKERRQ(ierr);
113778a2246SBarry Smith     ierr = VecCreateMPIWithArray(red->psubcomm->dupparent,1,mloc_sub,PETSC_DECIDE,PETSC_NULL,&red->ydup);CHKERRQ(ierr);
1143f457be1SHong Zhang 
1153f457be1SHong Zhang     /* create vec scatters */
1163f457be1SHong Zhang     if (!red->scatterin) {
1173f457be1SHong Zhang       IS       is1,is2;
1183f457be1SHong Zhang       PetscInt *idx1,*idx2,i,j,k;
11945fc02eaSBarry Smith 
1201d79065fSBarry Smith       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;
1253f457be1SHong 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);
1303f457be1SHong Zhang       ierr = VecScatterCreate(vec,is1,red->xdup,is2,&red->scatterin);CHKERRQ(ierr);
131fcfd50ebSBarry Smith       ierr = ISDestroy(&is1);CHKERRQ(ierr);
132fcfd50ebSBarry Smith       ierr = ISDestroy(&is2);CHKERRQ(ierr);
1333f457be1SHong Zhang 
1341fbd8f88SHong 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);
1363f457be1SHong Zhang       ierr = VecScatterCreate(red->xdup,is1,vec,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);
1404b9ad928SBarry Smith     }
1414b9ad928SBarry Smith   }
1426bf464f9SBarry Smith   ierr = VecDestroy(&vec);CHKERRQ(ierr);
1434b9ad928SBarry Smith 
1444b9ad928SBarry Smith   /* if pmatrix set by user is sequential then we do not need to gather the parallel matrix */
1453f457be1SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
146*2fa5cd67SKarl Rupp   if (size == 1) red->useparallelmat = PETSC_FALSE;
1474b9ad928SBarry Smith 
1484b9ad928SBarry Smith   if (red->useparallelmat) {
1494b9ad928SBarry Smith     if (pc->setupcalled == 1 && pc->flag == DIFFERENT_NONZERO_PATTERN) {
1504b9ad928SBarry Smith       /* destroy old matrices */
1516bf464f9SBarry Smith       ierr = MatDestroy(&red->pmats);CHKERRQ(ierr);
1524b9ad928SBarry Smith     } else if (pc->setupcalled == 1) {
1534b9ad928SBarry Smith       reuse = MAT_REUSE_MATRIX;
1544b9ad928SBarry Smith       str   = SAME_NONZERO_PATTERN;
1554b9ad928SBarry Smith     }
1564b9ad928SBarry Smith 
1573f457be1SHong Zhang     /* grab the parallel matrix and put it into processors of a subcomminicator */
158f664ae05SHong Zhang     /*--------------------------------------------------------------------------*/
159f664ae05SHong Zhang     ierr = VecGetLocalSize(red->ysub,&mlocal_sub);CHKERRQ(ierr);
16069db28dcSHong Zhang     ierr = MatGetRedundantMatrix(pc->pmat,red->psubcomm->n,red->psubcomm->comm,mlocal_sub,reuse,&red->pmats);CHKERRQ(ierr);
1613f457be1SHong Zhang     /* tell PC of the subcommunicator its operators */
16290f1c854SHong Zhang     ierr = KSPSetOperators(red->ksp,red->pmats,red->pmats,str);CHKERRQ(ierr);
1634b9ad928SBarry Smith   } else {
16490f1c854SHong Zhang     ierr = KSPSetOperators(red->ksp,pc->mat,pc->pmat,pc->flag);CHKERRQ(ierr);
1654b9ad928SBarry Smith   }
1660c24e6a1SHong Zhang   if (pc->setfromoptionscalled) {
1673e065800SHong Zhang     ierr = KSPSetFromOptions(red->ksp);CHKERRQ(ierr);
1680c24e6a1SHong Zhang   }
1693e065800SHong Zhang   ierr = KSPSetUp(red->ksp);CHKERRQ(ierr);
1704b9ad928SBarry Smith   PetscFunctionReturn(0);
1714b9ad928SBarry Smith }
1724b9ad928SBarry Smith 
1734b9ad928SBarry Smith #undef __FUNCT__
1744b9ad928SBarry Smith #define __FUNCT__ "PCApply_Redundant"
1756849ba73SBarry Smith static PetscErrorCode PCApply_Redundant(PC pc,Vec x,Vec y)
1764b9ad928SBarry Smith {
1774b9ad928SBarry Smith   PC_Redundant   *red = (PC_Redundant*)pc->data;
178dfbe8321SBarry Smith   PetscErrorCode ierr;
1793f457be1SHong Zhang   PetscScalar    *array;
1804b9ad928SBarry Smith 
1814b9ad928SBarry Smith   PetscFunctionBegin;
1823f457be1SHong Zhang   /* scatter x to xdup */
183ca9f406cSSatish Balay   ierr = VecScatterBegin(red->scatterin,x,red->xdup,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
184ca9f406cSSatish Balay   ierr = VecScatterEnd(red->scatterin,x,red->xdup,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1853f457be1SHong Zhang 
1863f457be1SHong Zhang   /* place xdup's local array into xsub */
1873f457be1SHong Zhang   ierr = VecGetArray(red->xdup,&array);CHKERRQ(ierr);
1883f457be1SHong Zhang   ierr = VecPlaceArray(red->xsub,(const PetscScalar*)array);CHKERRQ(ierr);
1894b9ad928SBarry Smith 
1904b9ad928SBarry Smith   /* apply preconditioner on each processor */
19183ab6a24SBarry Smith   ierr = KSPSolve(red->ksp,red->xsub,red->ysub);CHKERRQ(ierr);
1923f457be1SHong Zhang   ierr = VecResetArray(red->xsub);CHKERRQ(ierr);
1933f457be1SHong Zhang   ierr = VecRestoreArray(red->xdup,&array);CHKERRQ(ierr);
1944b9ad928SBarry Smith 
1953f457be1SHong Zhang   /* place ysub's local array into ydup */
1963f457be1SHong Zhang   ierr = VecGetArray(red->ysub,&array);CHKERRQ(ierr);
1973f457be1SHong Zhang   ierr = VecPlaceArray(red->ydup,(const PetscScalar*)array);CHKERRQ(ierr);
1983f457be1SHong Zhang 
1993f457be1SHong Zhang   /* scatter ydup to y */
200ca9f406cSSatish Balay   ierr = VecScatterBegin(red->scatterout,red->ydup,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
201ca9f406cSSatish Balay   ierr = VecScatterEnd(red->scatterout,red->ydup,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
2023f457be1SHong Zhang   ierr = VecResetArray(red->ydup);CHKERRQ(ierr);
2033f457be1SHong Zhang   ierr = VecRestoreArray(red->ysub,&array);CHKERRQ(ierr);
2044b9ad928SBarry Smith   PetscFunctionReturn(0);
2054b9ad928SBarry Smith }
2064b9ad928SBarry Smith 
2074b9ad928SBarry Smith #undef __FUNCT__
2081ea5a559SBarry Smith #define __FUNCT__ "PCReset_Redundant"
2091ea5a559SBarry Smith static PetscErrorCode PCReset_Redundant(PC pc)
2104b9ad928SBarry Smith {
2114b9ad928SBarry Smith   PC_Redundant   *red = (PC_Redundant*)pc->data;
212dfbe8321SBarry Smith   PetscErrorCode ierr;
2134b9ad928SBarry Smith 
2144b9ad928SBarry Smith   PetscFunctionBegin;
2156bf464f9SBarry Smith   ierr = VecScatterDestroy(&red->scatterin);CHKERRQ(ierr);
2166bf464f9SBarry Smith   ierr = VecScatterDestroy(&red->scatterout);CHKERRQ(ierr);
2176bf464f9SBarry Smith   ierr = VecDestroy(&red->ysub);CHKERRQ(ierr);
2186bf464f9SBarry Smith   ierr = VecDestroy(&red->xsub);CHKERRQ(ierr);
2196bf464f9SBarry Smith   ierr = VecDestroy(&red->xdup);CHKERRQ(ierr);
2206bf464f9SBarry Smith   ierr = VecDestroy(&red->ydup);CHKERRQ(ierr);
2216bf464f9SBarry Smith   ierr = MatDestroy(&red->pmats);CHKERRQ(ierr);
222546c30cbSLisandro Dalcin   if (red->ksp) {ierr = KSPReset(red->ksp);CHKERRQ(ierr);}
2231ea5a559SBarry Smith   PetscFunctionReturn(0);
2241ea5a559SBarry Smith }
2251ea5a559SBarry Smith 
2261ea5a559SBarry Smith #undef __FUNCT__
2271ea5a559SBarry Smith #define __FUNCT__ "PCDestroy_Redundant"
2281ea5a559SBarry Smith static PetscErrorCode PCDestroy_Redundant(PC pc)
2291ea5a559SBarry Smith {
2301ea5a559SBarry Smith   PC_Redundant   *red = (PC_Redundant*)pc->data;
2311ea5a559SBarry Smith   PetscErrorCode ierr;
2321ea5a559SBarry Smith 
2331ea5a559SBarry Smith   PetscFunctionBegin;
2341ea5a559SBarry Smith   ierr = PCReset_Redundant(pc);CHKERRQ(ierr);
2356bf464f9SBarry Smith   ierr = KSPDestroy(&red->ksp);CHKERRQ(ierr);
2366bf464f9SBarry Smith   ierr = PetscSubcommDestroy(&red->psubcomm);CHKERRQ(ierr);
237c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
2384b9ad928SBarry Smith   PetscFunctionReturn(0);
2394b9ad928SBarry Smith }
2404b9ad928SBarry Smith 
2414b9ad928SBarry Smith #undef __FUNCT__
2424b9ad928SBarry Smith #define __FUNCT__ "PCSetFromOptions_Redundant"
2436849ba73SBarry Smith static PetscErrorCode PCSetFromOptions_Redundant(PC pc)
2444b9ad928SBarry Smith {
245a98ce0f4SHong Zhang   PetscErrorCode ierr;
246a98ce0f4SHong Zhang   PC_Redundant   *red = (PC_Redundant*)pc->data;
247a98ce0f4SHong Zhang 
2484b9ad928SBarry Smith   PetscFunctionBegin;
249a98ce0f4SHong Zhang   ierr = PetscOptionsHead("Redundant options");CHKERRQ(ierr);
25009a6bc64SHong Zhang   ierr = PetscOptionsInt("-pc_redundant_number","Number of redundant pc","PCRedundantSetNumber",red->nsubcomm,&red->nsubcomm,0);CHKERRQ(ierr);
251a98ce0f4SHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
2524b9ad928SBarry Smith   PetscFunctionReturn(0);
2534b9ad928SBarry Smith }
2544b9ad928SBarry Smith 
2554b9ad928SBarry Smith EXTERN_C_BEGIN
2564b9ad928SBarry Smith #undef __FUNCT__
25709a6bc64SHong Zhang #define __FUNCT__ "PCRedundantSetNumber_Redundant"
2587087cfbeSBarry Smith 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 EXTERN_C_END
26709a6bc64SHong Zhang 
26809a6bc64SHong Zhang #undef __FUNCT__
26909a6bc64SHong Zhang #define __FUNCT__ "PCRedundantSetNumber"
27009a6bc64SHong Zhang /*@
27109a6bc64SHong Zhang    PCRedundantSetNumber - Sets the number of redundant preconditioner contexts.
27209a6bc64SHong Zhang 
2733f9fe445SBarry Smith    Logically Collective on PC
27409a6bc64SHong Zhang 
27509a6bc64SHong Zhang    Input Parameters:
27609a6bc64SHong Zhang +  pc - the preconditioner context
2779b21d695SBarry Smith -  nredundant - number of redundant preconditioner contexts; for example if you are using 64 MPI processes and
2789b21d695SBarry Smith                               use an nredundant of 4 there will be 4 parallel solves each on 16 = 64/4 processes.
27909a6bc64SHong Zhang 
28009a6bc64SHong Zhang    Level: advanced
28109a6bc64SHong Zhang 
28209a6bc64SHong Zhang .keywords: PC, redundant solve
28309a6bc64SHong Zhang @*/
2847087cfbeSBarry Smith PetscErrorCode  PCRedundantSetNumber(PC pc,PetscInt nredundant)
28509a6bc64SHong Zhang {
2864ac538c5SBarry Smith   PetscErrorCode ierr;
28709a6bc64SHong Zhang 
28809a6bc64SHong Zhang   PetscFunctionBegin;
2890700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
29065e19b50SBarry Smith   if (nredundant <= 0) SETERRQ1(((PetscObject)pc)->comm,PETSC_ERR_ARG_WRONG, "num of redundant pc %D must be positive",nredundant);
2914ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCRedundantSetNumber_C",(PC,PetscInt),(pc,nredundant));CHKERRQ(ierr);
29209a6bc64SHong Zhang   PetscFunctionReturn(0);
29309a6bc64SHong Zhang }
29409a6bc64SHong Zhang 
29509a6bc64SHong Zhang EXTERN_C_BEGIN
29609a6bc64SHong Zhang #undef __FUNCT__
2974b9ad928SBarry Smith #define __FUNCT__ "PCRedundantSetScatter_Redundant"
2987087cfbeSBarry Smith PetscErrorCode  PCRedundantSetScatter_Redundant(PC pc,VecScatter in,VecScatter out)
2994b9ad928SBarry Smith {
3004b9ad928SBarry Smith   PC_Redundant   *red = (PC_Redundant*)pc->data;
301dfbe8321SBarry Smith   PetscErrorCode ierr;
3024b9ad928SBarry Smith 
3034b9ad928SBarry Smith   PetscFunctionBegin;
3044b9ad928SBarry Smith   ierr = PetscObjectReference((PetscObject)in);CHKERRQ(ierr);
3056bf464f9SBarry Smith   ierr = VecScatterDestroy(&red->scatterin);CHKERRQ(ierr);
306*2fa5cd67SKarl Rupp 
307c3122656SLisandro Dalcin   red->scatterin  = in;
308*2fa5cd67SKarl Rupp 
3094b9ad928SBarry Smith   ierr = PetscObjectReference((PetscObject)out);CHKERRQ(ierr);
3106bf464f9SBarry Smith   ierr = VecScatterDestroy(&red->scatterout);CHKERRQ(ierr);
311c3122656SLisandro Dalcin   red->scatterout = out;
3124b9ad928SBarry Smith   PetscFunctionReturn(0);
3134b9ad928SBarry Smith }
3144b9ad928SBarry Smith EXTERN_C_END
3154b9ad928SBarry Smith 
3164b9ad928SBarry Smith #undef __FUNCT__
3174b9ad928SBarry Smith #define __FUNCT__ "PCRedundantSetScatter"
3184b9ad928SBarry Smith /*@
3194b9ad928SBarry Smith    PCRedundantSetScatter - Sets the scatter used to copy values into the
3204b9ad928SBarry Smith      redundant local solve and the scatter to move them back into the global
3214b9ad928SBarry Smith      vector.
3224b9ad928SBarry Smith 
3233f9fe445SBarry Smith    Logically Collective on PC
3244b9ad928SBarry Smith 
3254b9ad928SBarry Smith    Input Parameters:
3264b9ad928SBarry Smith +  pc - the preconditioner context
3274b9ad928SBarry Smith .  in - the scatter to move the values in
3284b9ad928SBarry Smith -  out - the scatter to move them out
3294b9ad928SBarry Smith 
3304b9ad928SBarry Smith    Level: advanced
3314b9ad928SBarry Smith 
3324b9ad928SBarry Smith .keywords: PC, redundant solve
3334b9ad928SBarry Smith @*/
3347087cfbeSBarry Smith PetscErrorCode  PCRedundantSetScatter(PC pc,VecScatter in,VecScatter out)
3354b9ad928SBarry Smith {
3364ac538c5SBarry Smith   PetscErrorCode ierr;
3374b9ad928SBarry Smith 
3384b9ad928SBarry Smith   PetscFunctionBegin;
3390700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
3400700a824SBarry Smith   PetscValidHeaderSpecific(in,VEC_SCATTER_CLASSID,2);
3410700a824SBarry Smith   PetscValidHeaderSpecific(out,VEC_SCATTER_CLASSID,3);
3424ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCRedundantSetScatter_C",(PC,VecScatter,VecScatter),(pc,in,out));CHKERRQ(ierr);
3434b9ad928SBarry Smith   PetscFunctionReturn(0);
3444b9ad928SBarry Smith }
3454b9ad928SBarry Smith 
3464b9ad928SBarry Smith EXTERN_C_BEGIN
3474b9ad928SBarry Smith #undef __FUNCT__
34883ab6a24SBarry Smith #define __FUNCT__ "PCRedundantGetKSP_Redundant"
34983ab6a24SBarry Smith 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;
366*2fa5cd67SKarl 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 EXTERN_C_END
3824b9ad928SBarry Smith 
3834b9ad928SBarry Smith #undef __FUNCT__
38483ab6a24SBarry Smith #define __FUNCT__ "PCRedundantGetKSP"
3854b9ad928SBarry Smith /*@
38683ab6a24SBarry Smith    PCRedundantGetKSP - Gets the less parallel KSP created by the redundant PC.
3874b9ad928SBarry Smith 
3884b9ad928SBarry Smith    Not Collective
3894b9ad928SBarry Smith 
3904b9ad928SBarry Smith    Input Parameter:
3914b9ad928SBarry Smith .  pc - the preconditioner context
3924b9ad928SBarry Smith 
3934b9ad928SBarry Smith    Output Parameter:
39483ab6a24SBarry Smith .  innerksp - the KSP on the smaller set of processes
3954b9ad928SBarry Smith 
3964b9ad928SBarry Smith    Level: advanced
3974b9ad928SBarry Smith 
3984b9ad928SBarry Smith .keywords: PC, redundant solve
3994b9ad928SBarry Smith @*/
40083ab6a24SBarry Smith PetscErrorCode  PCRedundantGetKSP(PC pc,KSP *innerksp)
4014b9ad928SBarry Smith {
4024ac538c5SBarry Smith   PetscErrorCode ierr;
4034b9ad928SBarry Smith 
4044b9ad928SBarry Smith   PetscFunctionBegin;
4050700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
40683ab6a24SBarry Smith   PetscValidPointer(innerksp,2);
40783ab6a24SBarry Smith   ierr = PetscTryMethod(pc,"PCRedundantGetKSP_C",(PC,KSP*),(pc,innerksp));CHKERRQ(ierr);
4084b9ad928SBarry Smith   PetscFunctionReturn(0);
4094b9ad928SBarry Smith }
4104b9ad928SBarry Smith 
4114b9ad928SBarry Smith EXTERN_C_BEGIN
4124b9ad928SBarry Smith #undef __FUNCT__
4134b9ad928SBarry Smith #define __FUNCT__ "PCRedundantGetOperators_Redundant"
4147087cfbeSBarry Smith PetscErrorCode  PCRedundantGetOperators_Redundant(PC pc,Mat *mat,Mat *pmat)
4154b9ad928SBarry Smith {
4164b9ad928SBarry Smith   PC_Redundant *red = (PC_Redundant*)pc->data;
4174b9ad928SBarry Smith 
4184b9ad928SBarry Smith   PetscFunctionBegin;
419b3804887SHong Zhang   if (mat)  *mat  = red->pmats;
420b3804887SHong Zhang   if (pmat) *pmat = red->pmats;
4214b9ad928SBarry Smith   PetscFunctionReturn(0);
4224b9ad928SBarry Smith }
4234b9ad928SBarry Smith EXTERN_C_END
4244b9ad928SBarry Smith 
4254b9ad928SBarry Smith #undef __FUNCT__
4264b9ad928SBarry Smith #define __FUNCT__ "PCRedundantGetOperators"
4274b9ad928SBarry Smith /*@
4284b9ad928SBarry Smith    PCRedundantGetOperators - gets the sequential matrix and preconditioner matrix
4294b9ad928SBarry Smith 
4304b9ad928SBarry Smith    Not Collective
4314b9ad928SBarry Smith 
4324b9ad928SBarry Smith    Input Parameter:
4334b9ad928SBarry Smith .  pc - the preconditioner context
4344b9ad928SBarry Smith 
4354b9ad928SBarry Smith    Output Parameters:
4364b9ad928SBarry Smith +  mat - the matrix
4374b9ad928SBarry Smith -  pmat - the (possibly different) preconditioner matrix
4384b9ad928SBarry Smith 
4394b9ad928SBarry Smith    Level: advanced
4404b9ad928SBarry Smith 
4414b9ad928SBarry Smith .keywords: PC, redundant solve
4424b9ad928SBarry Smith @*/
4437087cfbeSBarry Smith PetscErrorCode  PCRedundantGetOperators(PC pc,Mat *mat,Mat *pmat)
4444b9ad928SBarry Smith {
4454ac538c5SBarry Smith   PetscErrorCode ierr;
4464b9ad928SBarry Smith 
4474b9ad928SBarry Smith   PetscFunctionBegin;
4480700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
4494482741eSBarry Smith   if (mat)  PetscValidPointer(mat,2);
4504482741eSBarry Smith   if (pmat) PetscValidPointer(pmat,3);
4514ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCRedundantGetOperators_C",(PC,Mat*,Mat*),(pc,mat,pmat));CHKERRQ(ierr);
4524b9ad928SBarry Smith   PetscFunctionReturn(0);
4534b9ad928SBarry Smith }
4544b9ad928SBarry Smith 
4554b9ad928SBarry Smith /* -------------------------------------------------------------------------------------*/
45637a17b4dSBarry Smith /*MC
45783ab6a24SBarry Smith      PCREDUNDANT - Runs a KSP solver with preconditioner for the entire problem on subgroups of processors
45837a17b4dSBarry Smith 
45983ab6a24SBarry Smith      Options for the redundant preconditioners can be set with -redundant_pc_xxx for the redundant KSP with -redundant_ksp_xxx
46037a17b4dSBarry Smith 
46109391456SBarry Smith   Options Database:
4629b21d695SBarry Smith .  -pc_redundant_number <n> - number of redundant solves, for example if you are using 64 MPI processes and
4639b21d695SBarry Smith                               use an n of 4 there will be 4 parallel solves each on 16 = 64/4 processes.
46409391456SBarry Smith 
46537a17b4dSBarry Smith    Level: intermediate
46637a17b4dSBarry Smith 
46783ab6a24SBarry Smith    Notes: The default KSP is preonly and the default PC is LU.
46883ab6a24SBarry Smith 
46983ab6a24SBarry Smith    Developer Notes: Note that PCSetInitialGuessNonzero()  is not used by this class but likely should be.
4709cfaa89bSBarry Smith 
47137a17b4dSBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PCRedundantSetScatter(),
47283ab6a24SBarry Smith            PCRedundantGetKSP(), PCRedundantGetOperators(), PCRedundantSetNumber()
47337a17b4dSBarry Smith M*/
47437a17b4dSBarry Smith 
4754b9ad928SBarry Smith EXTERN_C_BEGIN
4764b9ad928SBarry Smith #undef __FUNCT__
4774b9ad928SBarry Smith #define __FUNCT__ "PCCreate_Redundant"
4787087cfbeSBarry Smith PetscErrorCode  PCCreate_Redundant(PC pc)
4794b9ad928SBarry Smith {
480dfbe8321SBarry Smith   PetscErrorCode ierr;
4814b9ad928SBarry Smith   PC_Redundant   *red;
48269db28dcSHong Zhang   PetscMPIInt    size;
4833f457be1SHong Zhang 
4844b9ad928SBarry Smith   PetscFunctionBegin;
48538f2d2fdSLisandro Dalcin   ierr = PetscNewLog(pc,PC_Redundant,&red);CHKERRQ(ierr);
4867adad957SLisandro Dalcin   ierr = MPI_Comm_size(((PetscObject)pc)->comm,&size);CHKERRQ(ierr);
487*2fa5cd67SKarl Rupp 
48869db28dcSHong Zhang   red->nsubcomm       = size;
4894b9ad928SBarry Smith   red->useparallelmat = PETSC_TRUE;
4901fbd8f88SHong Zhang   pc->data            = (void*)red;
4914b9ad928SBarry Smith 
4924b9ad928SBarry Smith   pc->ops->apply          = PCApply_Redundant;
4934b9ad928SBarry Smith   pc->ops->applytranspose = 0;
4944b9ad928SBarry Smith   pc->ops->setup          = PCSetUp_Redundant;
4954b9ad928SBarry Smith   pc->ops->destroy        = PCDestroy_Redundant;
4961ea5a559SBarry Smith   pc->ops->reset          = PCReset_Redundant;
4974b9ad928SBarry Smith   pc->ops->setfromoptions = PCSetFromOptions_Redundant;
4984b9ad928SBarry Smith   pc->ops->view           = PCView_Redundant;
499*2fa5cd67SKarl Rupp 
5004b9ad928SBarry Smith   ierr                    = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCRedundantSetScatter_C","PCRedundantSetScatter_Redundant",
5014b9ad928SBarry Smith                                                               PCRedundantSetScatter_Redundant);CHKERRQ(ierr);
50209a6bc64SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCRedundantSetNumber_C","PCRedundantSetNumber_Redundant",
50309a6bc64SHong Zhang                                            PCRedundantSetNumber_Redundant);CHKERRQ(ierr);
50483ab6a24SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCRedundantGetKSP_C","PCRedundantGetKSP_Redundant",
50583ab6a24SBarry Smith                                            PCRedundantGetKSP_Redundant);CHKERRQ(ierr);
5064b9ad928SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCRedundantGetOperators_C","PCRedundantGetOperators_Redundant",
5074b9ad928SBarry Smith                                            PCRedundantGetOperators_Redundant);CHKERRQ(ierr);
5084b9ad928SBarry Smith   PetscFunctionReturn(0);
5094b9ad928SBarry Smith }
5104b9ad928SBarry Smith EXTERN_C_END
511