xref: /petsc/src/ksp/pc/impls/redundant/redundant.c (revision e5acf8a457bbca6d488ec50ccac76a0142f931f5)
1dba47a55SKris Buschelman 
24b9ad928SBarry Smith /*
33f457be1SHong Zhang   This file defines a "solve the problem redundantly on each subgroup of processor" preconditioner.
44b9ad928SBarry Smith */
5af0996ceSBarry 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);
373f08860eSBarry Smith       ierr = PetscViewerGetSubViewer(viewer,((PetscObject)red->pc)->comm,&subviewer);CHKERRQ(ierr);
38f5dd71faSBarry Smith       if (!red->psubcomm->color) { /* only view first redundant pc */
391575c14dSBarry Smith         ierr = PetscViewerASCIIPushTab(subviewer);CHKERRQ(ierr);
403e065800SHong Zhang         ierr = KSPView(red->ksp,subviewer);CHKERRQ(ierr);
411575c14dSBarry Smith         ierr = PetscViewerASCIIPopTab(subviewer);CHKERRQ(ierr);
424b9ad928SBarry Smith       }
433f08860eSBarry Smith       ierr = PetscViewerRestoreSubViewer(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;
613f457be1SHong Zhang 
624b9ad928SBarry Smith   PetscFunctionBegin;
63ce94432eSBarry Smith   ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
64ddc54837SHong Zhang 
65ddc54837SHong Zhang   /* if pmatrix set by user is sequential then we do not need to gather the parallel matrix */
66ddc54837SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
67ddc54837SHong Zhang   if (size == 1) red->useparallelmat = PETSC_FALSE;
681fbd8f88SHong Zhang 
694b9ad928SBarry Smith   if (!pc->setupcalled) {
701b81debcSHong Zhang     PetscInt mloc_sub;
7175024027SHong Zhang     if (!red->psubcomm) { /* create red->psubcomm, new ksp and pc over subcomm */
7275024027SHong Zhang       KSP ksp;
7375024027SHong Zhang       ierr = PCRedundantGetKSP(pc,&ksp);CHKERRQ(ierr);
74306c2d5bSBarry Smith       subcomm = PetscSubcommChild(red->psubcomm);
751b81debcSHong Zhang     }
7675024027SHong Zhang     subcomm = PetscSubcommChild(red->psubcomm);
771fbd8f88SHong Zhang 
781b81debcSHong Zhang     if (red->useparallelmat) {
791b81debcSHong Zhang       /* grab the parallel matrix and put it into processors of a subcomminicator */
8053cd1579SHong Zhang       ierr = MatCreateRedundantMatrix(pc->pmat,red->psubcomm->n,subcomm,MAT_INITIAL_MATRIX,&red->pmats);CHKERRQ(ierr);
8123ee1639SBarry Smith       ierr = KSPSetOperators(red->ksp,red->pmats,red->pmats);CHKERRQ(ierr);
824b9ad928SBarry Smith 
831b81debcSHong Zhang       /* get working vectors xsub and ysub */
842a7a6963SBarry Smith       ierr = MatCreateVecs(red->pmats,&red->xsub,&red->ysub);CHKERRQ(ierr);
852fa5cd67SKarl Rupp 
868b96b0d2SHong Zhang       /* create working vectors xdup and ydup.
878b96b0d2SHong Zhang        xdup concatenates all xsub's contigously to form a mpi vector over dupcomm  (see PetscSubcommCreate_interlaced())
888b96b0d2SHong Zhang        ydup concatenates all ysub and has empty local arrays because ysub's arrays will be place into it.
89ce94432eSBarry Smith        Note: we use communicator dupcomm, not PetscObjectComm((PetscObject)pc)! */
901b81debcSHong Zhang       ierr = MatGetLocalSize(red->pmats,&mloc_sub,NULL);CHKERRQ(ierr);
9136be1a5eSBarry Smith       ierr = VecCreateMPI(PetscSubcommContiguousParent(red->psubcomm),mloc_sub,PETSC_DECIDE,&red->xdup);CHKERRQ(ierr);
9236be1a5eSBarry Smith       ierr = VecCreateMPIWithArray(PetscSubcommContiguousParent(red->psubcomm),1,mloc_sub,PETSC_DECIDE,NULL,&red->ydup);CHKERRQ(ierr);
933f457be1SHong Zhang 
94f68be91cSHong Zhang       /* create vecscatters */
95f68be91cSHong Zhang       if (!red->scatterin) { /* efficiency of scatterin is independent from psubcomm_type! */
963f457be1SHong Zhang         IS       is1,is2;
973f457be1SHong Zhang         PetscInt *idx1,*idx2,i,j,k;
9845fc02eaSBarry Smith 
992a7a6963SBarry Smith         ierr = MatCreateVecs(pc->pmat,&x,0);CHKERRQ(ierr);
1001b81debcSHong Zhang         ierr = VecGetSize(x,&M);CHKERRQ(ierr);
1011b81debcSHong Zhang         ierr = VecGetOwnershipRange(x,&mstart,&mend);CHKERRQ(ierr);
1021b81debcSHong Zhang         mlocal = mend - mstart;
103dcca6d9dSJed Brown         ierr = PetscMalloc2(red->psubcomm->n*mlocal,&idx1,red->psubcomm->n*mlocal,&idx2);CHKERRQ(ierr);
1043f457be1SHong Zhang         j    = 0;
1051fbd8f88SHong Zhang         for (k=0; k<red->psubcomm->n; k++) {
1063f457be1SHong Zhang           for (i=mstart; i<mend; i++) {
1073f457be1SHong Zhang             idx1[j]   = i;
108ddc54837SHong Zhang             idx2[j++] = i + M*k;
1093f457be1SHong Zhang           }
1103f457be1SHong Zhang         }
11170b3c8c7SBarry Smith         ierr = ISCreateGeneral(comm,red->psubcomm->n*mlocal,idx1,PETSC_COPY_VALUES,&is1);CHKERRQ(ierr);
11270b3c8c7SBarry Smith         ierr = ISCreateGeneral(comm,red->psubcomm->n*mlocal,idx2,PETSC_COPY_VALUES,&is2);CHKERRQ(ierr);
113ddc54837SHong Zhang         ierr = VecScatterCreate(x,is1,red->xdup,is2,&red->scatterin);CHKERRQ(ierr);
114fcfd50ebSBarry Smith         ierr = ISDestroy(&is1);CHKERRQ(ierr);
115fcfd50ebSBarry Smith         ierr = ISDestroy(&is2);CHKERRQ(ierr);
1163f457be1SHong Zhang 
1176909748bSHong Zhang         /* Impl below is good for PETSC_SUBCOMM_INTERLACED (no inter-process communication) and PETSC_SUBCOMM_CONTIGUOUS (communication within subcomm) */
118ddc54837SHong Zhang         ierr = ISCreateStride(comm,mlocal,mstart+ red->psubcomm->color*M,1,&is1);CHKERRQ(ierr);
1193f457be1SHong Zhang         ierr = ISCreateStride(comm,mlocal,mstart,1,&is2);CHKERRQ(ierr);
120ddc54837SHong Zhang         ierr = VecScatterCreate(red->xdup,is1,x,is2,&red->scatterout);CHKERRQ(ierr);
121fcfd50ebSBarry Smith         ierr = ISDestroy(&is1);CHKERRQ(ierr);
122fcfd50ebSBarry Smith         ierr = ISDestroy(&is2);CHKERRQ(ierr);
1231d79065fSBarry Smith         ierr = PetscFree2(idx1,idx2);CHKERRQ(ierr);
124ddc54837SHong Zhang         ierr = VecDestroy(&x);CHKERRQ(ierr);
1251b81debcSHong Zhang       }
126ab661555SHong Zhang     } else { /* !red->useparallelmat */
12723ee1639SBarry Smith       ierr = KSPSetOperators(red->ksp,pc->mat,pc->pmat);CHKERRQ(ierr);
1281b81debcSHong Zhang     }
129ab661555SHong Zhang   } else { /* pc->setupcalled */
1304b9ad928SBarry Smith     if (red->useparallelmat) {
131ab661555SHong Zhang       MatReuse       reuse;
1321b81debcSHong Zhang       /* grab the parallel matrix and put it into processors of a subcomminicator */
1331b81debcSHong Zhang       /*--------------------------------------------------------------------------*/
134ab661555SHong Zhang       if (pc->flag == DIFFERENT_NONZERO_PATTERN) {
1354b9ad928SBarry Smith         /* destroy old matrices */
1366bf464f9SBarry Smith         ierr  = MatDestroy(&red->pmats);CHKERRQ(ierr);
137ab661555SHong Zhang         reuse = MAT_INITIAL_MATRIX;
1384b9ad928SBarry Smith       } else {
139ab661555SHong Zhang         reuse = MAT_REUSE_MATRIX;
140ab661555SHong Zhang       }
141306c2d5bSBarry Smith       ierr = MatCreateRedundantMatrix(pc->pmat,red->psubcomm->n,PetscSubcommChild(red->psubcomm),reuse,&red->pmats);CHKERRQ(ierr);
14223ee1639SBarry Smith       ierr = KSPSetOperators(red->ksp,red->pmats,red->pmats);CHKERRQ(ierr);
143ab661555SHong Zhang     } else { /* !red->useparallelmat */
14423ee1639SBarry Smith       ierr = KSPSetOperators(red->ksp,pc->mat,pc->pmat);CHKERRQ(ierr);
1454b9ad928SBarry Smith     }
146ab661555SHong Zhang   }
1471b81debcSHong Zhang 
1480c24e6a1SHong Zhang   if (pc->setfromoptionscalled) {
1493e065800SHong Zhang     ierr = KSPSetFromOptions(red->ksp);CHKERRQ(ierr);
1500c24e6a1SHong Zhang   }
1513e065800SHong Zhang   ierr = KSPSetUp(red->ksp);CHKERRQ(ierr);
1524b9ad928SBarry Smith   PetscFunctionReturn(0);
1534b9ad928SBarry Smith }
1544b9ad928SBarry Smith 
1554b9ad928SBarry Smith #undef __FUNCT__
1564b9ad928SBarry Smith #define __FUNCT__ "PCApply_Redundant"
1576849ba73SBarry Smith static PetscErrorCode PCApply_Redundant(PC pc,Vec x,Vec y)
1584b9ad928SBarry Smith {
1594b9ad928SBarry Smith   PC_Redundant   *red = (PC_Redundant*)pc->data;
160dfbe8321SBarry Smith   PetscErrorCode ierr;
1613f457be1SHong Zhang   PetscScalar    *array;
1624b9ad928SBarry Smith 
1634b9ad928SBarry Smith   PetscFunctionBegin;
164ddc54837SHong Zhang   if (!red->useparallelmat) {
165ddc54837SHong Zhang     ierr = KSPSolve(red->ksp,x,y);CHKERRQ(ierr);
166ddc54837SHong Zhang     PetscFunctionReturn(0);
167ddc54837SHong Zhang   }
168ddc54837SHong Zhang 
1693f457be1SHong Zhang   /* scatter x to xdup */
170ca9f406cSSatish Balay   ierr = VecScatterBegin(red->scatterin,x,red->xdup,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
171ca9f406cSSatish Balay   ierr = VecScatterEnd(red->scatterin,x,red->xdup,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1723f457be1SHong Zhang 
1733f457be1SHong Zhang   /* place xdup's local array into xsub */
1743f457be1SHong Zhang   ierr = VecGetArray(red->xdup,&array);CHKERRQ(ierr);
1753f457be1SHong Zhang   ierr = VecPlaceArray(red->xsub,(const PetscScalar*)array);CHKERRQ(ierr);
1764b9ad928SBarry Smith 
1774b9ad928SBarry Smith   /* apply preconditioner on each processor */
17883ab6a24SBarry Smith   ierr = KSPSolve(red->ksp,red->xsub,red->ysub);CHKERRQ(ierr);
1793f457be1SHong Zhang   ierr = VecResetArray(red->xsub);CHKERRQ(ierr);
1803f457be1SHong Zhang   ierr = VecRestoreArray(red->xdup,&array);CHKERRQ(ierr);
1814b9ad928SBarry Smith 
1823f457be1SHong Zhang   /* place ysub's local array into ydup */
1833f457be1SHong Zhang   ierr = VecGetArray(red->ysub,&array);CHKERRQ(ierr);
1843f457be1SHong Zhang   ierr = VecPlaceArray(red->ydup,(const PetscScalar*)array);CHKERRQ(ierr);
1853f457be1SHong Zhang 
1863f457be1SHong Zhang   /* scatter ydup to y */
187ca9f406cSSatish Balay   ierr = VecScatterBegin(red->scatterout,red->ydup,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
188ca9f406cSSatish Balay   ierr = VecScatterEnd(red->scatterout,red->ydup,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1893f457be1SHong Zhang   ierr = VecResetArray(red->ydup);CHKERRQ(ierr);
1903f457be1SHong Zhang   ierr = VecRestoreArray(red->ysub,&array);CHKERRQ(ierr);
1914b9ad928SBarry Smith   PetscFunctionReturn(0);
1924b9ad928SBarry Smith }
1934b9ad928SBarry Smith 
1944b9ad928SBarry Smith #undef __FUNCT__
195d88bfacbSStefano Zampini #define __FUNCT__ "PCApplyTranspose_Redundant"
196d88bfacbSStefano Zampini static PetscErrorCode PCApplyTranspose_Redundant(PC pc,Vec x,Vec y)
197d88bfacbSStefano Zampini {
198d88bfacbSStefano Zampini   PC_Redundant   *red = (PC_Redundant*)pc->data;
199d88bfacbSStefano Zampini   PetscErrorCode ierr;
200d88bfacbSStefano Zampini   PetscScalar    *array;
201d88bfacbSStefano Zampini 
202d88bfacbSStefano Zampini   PetscFunctionBegin;
203d88bfacbSStefano Zampini   if (!red->useparallelmat) {
204d88bfacbSStefano Zampini     ierr = KSPSolveTranspose(red->ksp,x,y);CHKERRQ(ierr);
205d88bfacbSStefano Zampini     PetscFunctionReturn(0);
206d88bfacbSStefano Zampini   }
207d88bfacbSStefano Zampini 
208d88bfacbSStefano Zampini   /* scatter x to xdup */
209d88bfacbSStefano Zampini   ierr = VecScatterBegin(red->scatterin,x,red->xdup,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
210d88bfacbSStefano Zampini   ierr = VecScatterEnd(red->scatterin,x,red->xdup,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
211d88bfacbSStefano Zampini 
212d88bfacbSStefano Zampini   /* place xdup's local array into xsub */
213d88bfacbSStefano Zampini   ierr = VecGetArray(red->xdup,&array);CHKERRQ(ierr);
214d88bfacbSStefano Zampini   ierr = VecPlaceArray(red->xsub,(const PetscScalar*)array);CHKERRQ(ierr);
215d88bfacbSStefano Zampini 
216d88bfacbSStefano Zampini   /* apply preconditioner on each processor */
217d88bfacbSStefano Zampini   ierr = KSPSolveTranspose(red->ksp,red->xsub,red->ysub);CHKERRQ(ierr);
218d88bfacbSStefano Zampini   ierr = VecResetArray(red->xsub);CHKERRQ(ierr);
219d88bfacbSStefano Zampini   ierr = VecRestoreArray(red->xdup,&array);CHKERRQ(ierr);
220d88bfacbSStefano Zampini 
221d88bfacbSStefano Zampini   /* place ysub's local array into ydup */
222d88bfacbSStefano Zampini   ierr = VecGetArray(red->ysub,&array);CHKERRQ(ierr);
223d88bfacbSStefano Zampini   ierr = VecPlaceArray(red->ydup,(const PetscScalar*)array);CHKERRQ(ierr);
224d88bfacbSStefano Zampini 
225d88bfacbSStefano Zampini   /* scatter ydup to y */
226d88bfacbSStefano Zampini   ierr = VecScatterBegin(red->scatterout,red->ydup,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
227d88bfacbSStefano Zampini   ierr = VecScatterEnd(red->scatterout,red->ydup,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
228d88bfacbSStefano Zampini   ierr = VecResetArray(red->ydup);CHKERRQ(ierr);
229d88bfacbSStefano Zampini   ierr = VecRestoreArray(red->ysub,&array);CHKERRQ(ierr);
230d88bfacbSStefano Zampini   PetscFunctionReturn(0);
231d88bfacbSStefano Zampini }
232d88bfacbSStefano Zampini 
233d88bfacbSStefano Zampini #undef __FUNCT__
2341ea5a559SBarry Smith #define __FUNCT__ "PCReset_Redundant"
2351ea5a559SBarry Smith static PetscErrorCode PCReset_Redundant(PC pc)
2364b9ad928SBarry Smith {
2374b9ad928SBarry Smith   PC_Redundant   *red = (PC_Redundant*)pc->data;
238dfbe8321SBarry Smith   PetscErrorCode ierr;
2394b9ad928SBarry Smith 
2404b9ad928SBarry Smith   PetscFunctionBegin;
2411b81debcSHong Zhang   if (red->useparallelmat) {
2426bf464f9SBarry Smith     ierr = VecScatterDestroy(&red->scatterin);CHKERRQ(ierr);
2436bf464f9SBarry Smith     ierr = VecScatterDestroy(&red->scatterout);CHKERRQ(ierr);
2446bf464f9SBarry Smith     ierr = VecDestroy(&red->ysub);CHKERRQ(ierr);
2456bf464f9SBarry Smith     ierr = VecDestroy(&red->xsub);CHKERRQ(ierr);
2466bf464f9SBarry Smith     ierr = VecDestroy(&red->xdup);CHKERRQ(ierr);
2476bf464f9SBarry Smith     ierr = VecDestroy(&red->ydup);CHKERRQ(ierr);
2481b81debcSHong Zhang   }
2496bf464f9SBarry Smith   ierr = MatDestroy(&red->pmats);CHKERRQ(ierr);
2501b81debcSHong Zhang   ierr = KSPReset(red->ksp);CHKERRQ(ierr);
2511ea5a559SBarry Smith   PetscFunctionReturn(0);
2521ea5a559SBarry Smith }
2531ea5a559SBarry Smith 
2541ea5a559SBarry Smith #undef __FUNCT__
2551ea5a559SBarry Smith #define __FUNCT__ "PCDestroy_Redundant"
2561ea5a559SBarry Smith static PetscErrorCode PCDestroy_Redundant(PC pc)
2571ea5a559SBarry Smith {
2581ea5a559SBarry Smith   PC_Redundant   *red = (PC_Redundant*)pc->data;
2591ea5a559SBarry Smith   PetscErrorCode ierr;
2601ea5a559SBarry Smith 
2611ea5a559SBarry Smith   PetscFunctionBegin;
2621ea5a559SBarry Smith   ierr = PCReset_Redundant(pc);CHKERRQ(ierr);
2636bf464f9SBarry Smith   ierr = KSPDestroy(&red->ksp);CHKERRQ(ierr);
2646bf464f9SBarry Smith   ierr = PetscSubcommDestroy(&red->psubcomm);CHKERRQ(ierr);
265c31cb41cSBarry Smith   ierr = PetscFree(pc->data);CHKERRQ(ierr);
2664b9ad928SBarry Smith   PetscFunctionReturn(0);
2674b9ad928SBarry Smith }
2684b9ad928SBarry Smith 
2694b9ad928SBarry Smith #undef __FUNCT__
2704b9ad928SBarry Smith #define __FUNCT__ "PCSetFromOptions_Redundant"
2714416b707SBarry Smith static PetscErrorCode PCSetFromOptions_Redundant(PetscOptionItems *PetscOptionsObject,PC pc)
2724b9ad928SBarry Smith {
273a98ce0f4SHong Zhang   PetscErrorCode ierr;
274a98ce0f4SHong Zhang   PC_Redundant   *red = (PC_Redundant*)pc->data;
275a98ce0f4SHong Zhang 
2764b9ad928SBarry Smith   PetscFunctionBegin;
277e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"Redundant options");CHKERRQ(ierr);
27809a6bc64SHong Zhang   ierr = PetscOptionsInt("-pc_redundant_number","Number of redundant pc","PCRedundantSetNumber",red->nsubcomm,&red->nsubcomm,0);CHKERRQ(ierr);
279a98ce0f4SHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
2804b9ad928SBarry Smith   PetscFunctionReturn(0);
2814b9ad928SBarry Smith }
2824b9ad928SBarry Smith 
2834b9ad928SBarry Smith #undef __FUNCT__
28409a6bc64SHong Zhang #define __FUNCT__ "PCRedundantSetNumber_Redundant"
285f7a08781SBarry Smith static PetscErrorCode PCRedundantSetNumber_Redundant(PC pc,PetscInt nreds)
28609a6bc64SHong Zhang {
28709a6bc64SHong Zhang   PC_Redundant *red = (PC_Redundant*)pc->data;
28809a6bc64SHong Zhang 
28909a6bc64SHong Zhang   PetscFunctionBegin;
29009a6bc64SHong Zhang   red->nsubcomm = nreds;
29109a6bc64SHong Zhang   PetscFunctionReturn(0);
29209a6bc64SHong Zhang }
29309a6bc64SHong Zhang 
29409a6bc64SHong Zhang #undef __FUNCT__
29509a6bc64SHong Zhang #define __FUNCT__ "PCRedundantSetNumber"
29609a6bc64SHong Zhang /*@
29709a6bc64SHong Zhang    PCRedundantSetNumber - Sets the number of redundant preconditioner contexts.
29809a6bc64SHong Zhang 
2993f9fe445SBarry Smith    Logically Collective on PC
30009a6bc64SHong Zhang 
30109a6bc64SHong Zhang    Input Parameters:
30209a6bc64SHong Zhang +  pc - the preconditioner context
3039b21d695SBarry Smith -  nredundant - number of redundant preconditioner contexts; for example if you are using 64 MPI processes and
3049b21d695SBarry Smith                               use an nredundant of 4 there will be 4 parallel solves each on 16 = 64/4 processes.
30509a6bc64SHong Zhang 
30609a6bc64SHong Zhang    Level: advanced
30709a6bc64SHong Zhang 
30809a6bc64SHong Zhang .keywords: PC, redundant solve
30909a6bc64SHong Zhang @*/
3107087cfbeSBarry Smith PetscErrorCode PCRedundantSetNumber(PC pc,PetscInt nredundant)
31109a6bc64SHong Zhang {
3124ac538c5SBarry Smith   PetscErrorCode ierr;
31309a6bc64SHong Zhang 
31409a6bc64SHong Zhang   PetscFunctionBegin;
3150700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
316ce94432eSBarry Smith   if (nredundant <= 0) SETERRQ1(PetscObjectComm((PetscObject)pc),PETSC_ERR_ARG_WRONG, "num of redundant pc %D must be positive",nredundant);
3174ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCRedundantSetNumber_C",(PC,PetscInt),(pc,nredundant));CHKERRQ(ierr);
31809a6bc64SHong Zhang   PetscFunctionReturn(0);
31909a6bc64SHong Zhang }
32009a6bc64SHong Zhang 
32109a6bc64SHong Zhang #undef __FUNCT__
3224b9ad928SBarry Smith #define __FUNCT__ "PCRedundantSetScatter_Redundant"
323f7a08781SBarry Smith static PetscErrorCode PCRedundantSetScatter_Redundant(PC pc,VecScatter in,VecScatter out)
3244b9ad928SBarry Smith {
3254b9ad928SBarry Smith   PC_Redundant   *red = (PC_Redundant*)pc->data;
326dfbe8321SBarry Smith   PetscErrorCode ierr;
3274b9ad928SBarry Smith 
3284b9ad928SBarry Smith   PetscFunctionBegin;
3294b9ad928SBarry Smith   ierr = PetscObjectReference((PetscObject)in);CHKERRQ(ierr);
3306bf464f9SBarry Smith   ierr = VecScatterDestroy(&red->scatterin);CHKERRQ(ierr);
3312fa5cd67SKarl Rupp 
332c3122656SLisandro Dalcin   red->scatterin  = in;
3332fa5cd67SKarl Rupp 
3344b9ad928SBarry Smith   ierr = PetscObjectReference((PetscObject)out);CHKERRQ(ierr);
3356bf464f9SBarry Smith   ierr = VecScatterDestroy(&red->scatterout);CHKERRQ(ierr);
336c3122656SLisandro Dalcin   red->scatterout = out;
3374b9ad928SBarry Smith   PetscFunctionReturn(0);
3384b9ad928SBarry Smith }
3394b9ad928SBarry Smith 
3404b9ad928SBarry Smith #undef __FUNCT__
3414b9ad928SBarry Smith #define __FUNCT__ "PCRedundantSetScatter"
3424b9ad928SBarry Smith /*@
3434b9ad928SBarry Smith    PCRedundantSetScatter - Sets the scatter used to copy values into the
3444b9ad928SBarry Smith      redundant local solve and the scatter to move them back into the global
3454b9ad928SBarry Smith      vector.
3464b9ad928SBarry Smith 
3473f9fe445SBarry Smith    Logically Collective on PC
3484b9ad928SBarry Smith 
3494b9ad928SBarry Smith    Input Parameters:
3504b9ad928SBarry Smith +  pc - the preconditioner context
3514b9ad928SBarry Smith .  in - the scatter to move the values in
3524b9ad928SBarry Smith -  out - the scatter to move them out
3534b9ad928SBarry Smith 
3544b9ad928SBarry Smith    Level: advanced
3554b9ad928SBarry Smith 
3564b9ad928SBarry Smith .keywords: PC, redundant solve
3574b9ad928SBarry Smith @*/
3587087cfbeSBarry Smith PetscErrorCode PCRedundantSetScatter(PC pc,VecScatter in,VecScatter out)
3594b9ad928SBarry Smith {
3604ac538c5SBarry Smith   PetscErrorCode ierr;
3614b9ad928SBarry Smith 
3624b9ad928SBarry Smith   PetscFunctionBegin;
3630700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
3640700a824SBarry Smith   PetscValidHeaderSpecific(in,VEC_SCATTER_CLASSID,2);
3650700a824SBarry Smith   PetscValidHeaderSpecific(out,VEC_SCATTER_CLASSID,3);
3664ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCRedundantSetScatter_C",(PC,VecScatter,VecScatter),(pc,in,out));CHKERRQ(ierr);
3674b9ad928SBarry Smith   PetscFunctionReturn(0);
3684b9ad928SBarry Smith }
3694b9ad928SBarry Smith 
3704b9ad928SBarry Smith #undef __FUNCT__
37183ab6a24SBarry Smith #define __FUNCT__ "PCRedundantGetKSP_Redundant"
372f7a08781SBarry Smith static PetscErrorCode PCRedundantGetKSP_Redundant(PC pc,KSP *innerksp)
3734b9ad928SBarry Smith {
3745f06b7aaSBarry Smith   PetscErrorCode ierr;
3754b9ad928SBarry Smith   PC_Redundant   *red = (PC_Redundant*)pc->data;
37675024027SHong Zhang   MPI_Comm       comm,subcomm;
37775024027SHong Zhang   const char     *prefix;
3784b9ad928SBarry Smith 
3794b9ad928SBarry Smith   PetscFunctionBegin;
38075024027SHong Zhang   if (!red->psubcomm) {
381*e5acf8a4SHong Zhang     ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
382*e5acf8a4SHong Zhang 
38375024027SHong Zhang     ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
38475024027SHong Zhang     ierr = PetscSubcommCreate(comm,&red->psubcomm);CHKERRQ(ierr);
38575024027SHong Zhang     ierr = PetscSubcommSetNumber(red->psubcomm,red->nsubcomm);CHKERRQ(ierr);
38675024027SHong Zhang     ierr = PetscSubcommSetType(red->psubcomm,PETSC_SUBCOMM_CONTIGUOUS);CHKERRQ(ierr);
387*e5acf8a4SHong Zhang 
388*e5acf8a4SHong Zhang     ierr = PetscSubcommSetOptionsPrefix(red->psubcomm,prefix);CHKERRQ(ierr);
389*e5acf8a4SHong Zhang     ierr = PetscSubcommSetFromOptions(red->psubcomm);CHKERRQ(ierr);
39075024027SHong Zhang     ierr = PetscLogObjectMemory((PetscObject)pc,sizeof(PetscSubcomm));CHKERRQ(ierr);
39175024027SHong Zhang 
39275024027SHong Zhang     /* create a new PC that processors in each subcomm have copy of */
39375024027SHong Zhang     subcomm = PetscSubcommChild(red->psubcomm);
39475024027SHong Zhang 
39575024027SHong Zhang     ierr = KSPCreate(subcomm,&red->ksp);CHKERRQ(ierr);
39675024027SHong Zhang     ierr = KSPSetErrorIfNotConverged(red->ksp,pc->erroriffailure);CHKERRQ(ierr);
39775024027SHong Zhang     ierr = PetscObjectIncrementTabLevel((PetscObject)red->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
39875024027SHong Zhang     ierr = PetscLogObjectParent((PetscObject)pc,(PetscObject)red->ksp);CHKERRQ(ierr);
39975024027SHong Zhang     ierr = KSPSetType(red->ksp,KSPPREONLY);CHKERRQ(ierr);
40075024027SHong Zhang     ierr = KSPGetPC(red->ksp,&red->pc);CHKERRQ(ierr);
40175024027SHong Zhang     ierr = PCSetType(red->pc,PCLU);CHKERRQ(ierr);
40275024027SHong Zhang 
40375024027SHong Zhang     ierr = KSPSetOptionsPrefix(red->ksp,prefix);CHKERRQ(ierr);
40475024027SHong Zhang     ierr = KSPAppendOptionsPrefix(red->ksp,"redundant_");CHKERRQ(ierr);
40575024027SHong Zhang   }
40683ab6a24SBarry Smith   *innerksp = red->ksp;
4074b9ad928SBarry Smith   PetscFunctionReturn(0);
4084b9ad928SBarry Smith }
4094b9ad928SBarry Smith 
4104b9ad928SBarry Smith #undef __FUNCT__
41183ab6a24SBarry Smith #define __FUNCT__ "PCRedundantGetKSP"
4124b9ad928SBarry Smith /*@
41383ab6a24SBarry Smith    PCRedundantGetKSP - Gets the less parallel KSP created by the redundant PC.
4144b9ad928SBarry Smith 
4154b9ad928SBarry Smith    Not Collective
4164b9ad928SBarry Smith 
4174b9ad928SBarry Smith    Input Parameter:
4184b9ad928SBarry Smith .  pc - the preconditioner context
4194b9ad928SBarry Smith 
4204b9ad928SBarry Smith    Output Parameter:
42183ab6a24SBarry Smith .  innerksp - the KSP on the smaller set of processes
4224b9ad928SBarry Smith 
4234b9ad928SBarry Smith    Level: advanced
4244b9ad928SBarry Smith 
4254b9ad928SBarry Smith .keywords: PC, redundant solve
4264b9ad928SBarry Smith @*/
42783ab6a24SBarry Smith PetscErrorCode PCRedundantGetKSP(PC pc,KSP *innerksp)
4284b9ad928SBarry Smith {
4294ac538c5SBarry Smith   PetscErrorCode ierr;
4304b9ad928SBarry Smith 
4314b9ad928SBarry Smith   PetscFunctionBegin;
4320700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
43383ab6a24SBarry Smith   PetscValidPointer(innerksp,2);
43483ab6a24SBarry Smith   ierr = PetscTryMethod(pc,"PCRedundantGetKSP_C",(PC,KSP*),(pc,innerksp));CHKERRQ(ierr);
4354b9ad928SBarry Smith   PetscFunctionReturn(0);
4364b9ad928SBarry Smith }
4374b9ad928SBarry Smith 
4384b9ad928SBarry Smith #undef __FUNCT__
4394b9ad928SBarry Smith #define __FUNCT__ "PCRedundantGetOperators_Redundant"
440f7a08781SBarry Smith static PetscErrorCode PCRedundantGetOperators_Redundant(PC pc,Mat *mat,Mat *pmat)
4414b9ad928SBarry Smith {
4424b9ad928SBarry Smith   PC_Redundant *red = (PC_Redundant*)pc->data;
4434b9ad928SBarry Smith 
4444b9ad928SBarry Smith   PetscFunctionBegin;
445b3804887SHong Zhang   if (mat)  *mat  = red->pmats;
446b3804887SHong Zhang   if (pmat) *pmat = red->pmats;
4474b9ad928SBarry Smith   PetscFunctionReturn(0);
4484b9ad928SBarry Smith }
4494b9ad928SBarry Smith 
4504b9ad928SBarry Smith #undef __FUNCT__
4514b9ad928SBarry Smith #define __FUNCT__ "PCRedundantGetOperators"
4524b9ad928SBarry Smith /*@
4534b9ad928SBarry Smith    PCRedundantGetOperators - gets the sequential matrix and preconditioner matrix
4544b9ad928SBarry Smith 
4554b9ad928SBarry Smith    Not Collective
4564b9ad928SBarry Smith 
4574b9ad928SBarry Smith    Input Parameter:
4584b9ad928SBarry Smith .  pc - the preconditioner context
4594b9ad928SBarry Smith 
4604b9ad928SBarry Smith    Output Parameters:
4614b9ad928SBarry Smith +  mat - the matrix
4624b9ad928SBarry Smith -  pmat - the (possibly different) preconditioner matrix
4634b9ad928SBarry Smith 
4644b9ad928SBarry Smith    Level: advanced
4654b9ad928SBarry Smith 
4664b9ad928SBarry Smith .keywords: PC, redundant solve
4674b9ad928SBarry Smith @*/
4687087cfbeSBarry Smith PetscErrorCode PCRedundantGetOperators(PC pc,Mat *mat,Mat *pmat)
4694b9ad928SBarry Smith {
4704ac538c5SBarry Smith   PetscErrorCode ierr;
4714b9ad928SBarry Smith 
4724b9ad928SBarry Smith   PetscFunctionBegin;
4730700a824SBarry Smith   PetscValidHeaderSpecific(pc,PC_CLASSID,1);
4744482741eSBarry Smith   if (mat)  PetscValidPointer(mat,2);
4754482741eSBarry Smith   if (pmat) PetscValidPointer(pmat,3);
4764ac538c5SBarry Smith   ierr = PetscTryMethod(pc,"PCRedundantGetOperators_C",(PC,Mat*,Mat*),(pc,mat,pmat));CHKERRQ(ierr);
4774b9ad928SBarry Smith   PetscFunctionReturn(0);
4784b9ad928SBarry Smith }
4794b9ad928SBarry Smith 
4804b9ad928SBarry Smith /* -------------------------------------------------------------------------------------*/
48137a17b4dSBarry Smith /*MC
48283ab6a24SBarry Smith      PCREDUNDANT - Runs a KSP solver with preconditioner for the entire problem on subgroups of processors
48337a17b4dSBarry Smith 
48483ab6a24SBarry Smith      Options for the redundant preconditioners can be set with -redundant_pc_xxx for the redundant KSP with -redundant_ksp_xxx
48537a17b4dSBarry Smith 
48609391456SBarry Smith   Options Database:
4879b21d695SBarry Smith .  -pc_redundant_number <n> - number of redundant solves, for example if you are using 64 MPI processes and
4889b21d695SBarry Smith                               use an n of 4 there will be 4 parallel solves each on 16 = 64/4 processes.
48909391456SBarry Smith 
49037a17b4dSBarry Smith    Level: intermediate
49137a17b4dSBarry Smith 
49283ab6a24SBarry Smith    Notes: The default KSP is preonly and the default PC is LU.
49383ab6a24SBarry Smith 
49483ab6a24SBarry Smith    Developer Notes: Note that PCSetInitialGuessNonzero()  is not used by this class but likely should be.
4959cfaa89bSBarry Smith 
49637a17b4dSBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PCRedundantSetScatter(),
49783ab6a24SBarry Smith            PCRedundantGetKSP(), PCRedundantGetOperators(), PCRedundantSetNumber()
49837a17b4dSBarry Smith M*/
49937a17b4dSBarry Smith 
5004b9ad928SBarry Smith #undef __FUNCT__
5014b9ad928SBarry Smith #define __FUNCT__ "PCCreate_Redundant"
5028cc058d9SJed Brown PETSC_EXTERN PetscErrorCode PCCreate_Redundant(PC pc)
5034b9ad928SBarry Smith {
504dfbe8321SBarry Smith   PetscErrorCode ierr;
5054b9ad928SBarry Smith   PC_Redundant   *red;
50669db28dcSHong Zhang   PetscMPIInt    size;
5073f457be1SHong Zhang 
5084b9ad928SBarry Smith   PetscFunctionBegin;
509b00a9115SJed Brown   ierr = PetscNewLog(pc,&red);CHKERRQ(ierr);
510ce94432eSBarry Smith   ierr = MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);CHKERRQ(ierr);
5112fa5cd67SKarl Rupp 
51269db28dcSHong Zhang   red->nsubcomm       = size;
5134b9ad928SBarry Smith   red->useparallelmat = PETSC_TRUE;
5141fbd8f88SHong Zhang   pc->data            = (void*)red;
5154b9ad928SBarry Smith 
5164b9ad928SBarry Smith   pc->ops->apply          = PCApply_Redundant;
517d88bfacbSStefano Zampini   pc->ops->applytranspose = PCApplyTranspose_Redundant;
5184b9ad928SBarry Smith   pc->ops->setup          = PCSetUp_Redundant;
5194b9ad928SBarry Smith   pc->ops->destroy        = PCDestroy_Redundant;
5201ea5a559SBarry Smith   pc->ops->reset          = PCReset_Redundant;
5214b9ad928SBarry Smith   pc->ops->setfromoptions = PCSetFromOptions_Redundant;
5224b9ad928SBarry Smith   pc->ops->view           = PCView_Redundant;
5232fa5cd67SKarl Rupp 
524bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCRedundantSetScatter_C",PCRedundantSetScatter_Redundant);CHKERRQ(ierr);
525bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCRedundantSetNumber_C",PCRedundantSetNumber_Redundant);CHKERRQ(ierr);
526bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCRedundantGetKSP_C",PCRedundantGetKSP_Redundant);CHKERRQ(ierr);
527bdf89e91SBarry Smith   ierr = PetscObjectComposeFunction((PetscObject)pc,"PCRedundantGetOperators_C",PCRedundantGetOperators_Redundant);CHKERRQ(ierr);
5284b9ad928SBarry Smith   PetscFunctionReturn(0);
5294b9ad928SBarry Smith }
530b2573a8aSBarry Smith 
531