xref: /petsc/src/ksp/pc/impls/redundant/redundant.c (revision 0c24e6a1e7d892dc3a9be60cfbbff2f1d2c50e14)
1dba47a55SKris Buschelman #define PETSCKSP_DLL
2dba47a55SKris Buschelman 
34b9ad928SBarry Smith /*
43f457be1SHong Zhang   This file defines a "solve the problem redundantly on each subgroup of processor" preconditioner.
54b9ad928SBarry Smith */
66356e834SBarry Smith #include "private/pcimpl.h"     /*I "petscpc.h" I*/
74b9ad928SBarry Smith #include "petscksp.h"
84b9ad928SBarry Smith 
94b9ad928SBarry Smith typedef struct {
104b9ad928SBarry Smith   PC           pc;                   /* actual preconditioner used on each processor */
113f457be1SHong Zhang   Vec          xsub,ysub;            /* vectors of a subcommunicator to hold parallel vectors of 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) */
154b9ad928SBarry Smith   PetscTruth   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;
2613f74950SBarry Smith   PetscMPIInt    rank;
2732077d6dSBarry Smith   PetscTruth     iascii,isstring;
28a47c9f9aSHong Zhang   PetscViewer    sviewer,subviewer;
291fbd8f88SHong Zhang   PetscInt       color = red->psubcomm->color;
304b9ad928SBarry Smith 
314b9ad928SBarry Smith   PetscFunctionBegin;
324b9ad928SBarry Smith   ierr = MPI_Comm_rank(pc->comm,&rank);CHKERRQ(ierr);
3332077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
344b9ad928SBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr);
3532077d6dSBarry Smith   if (iascii) {
36*0c24e6a1SHong Zhang     ierr = PetscViewerASCIIPrintf(viewer,"  Redundant solver preconditioner: First (color=0) of %D PCs follows\n",red->nsubcomm);CHKERRQ(ierr);
37a98ce0f4SHong Zhang     ierr = PetscViewerGetSubcomm(viewer,red->pc->comm,&subviewer);CHKERRQ(ierr);
38a98ce0f4SHong Zhang     if (!color) { /* only view first redundant pc */
394b9ad928SBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
40a47c9f9aSHong Zhang       ierr = PCView(red->pc,subviewer);CHKERRQ(ierr);
414b9ad928SBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
424b9ad928SBarry Smith     }
43a98ce0f4SHong Zhang     ierr = PetscViewerRestoreSubcomm(viewer,red->pc->comm,&subviewer);CHKERRQ(ierr);
44a98ce0f4SHong Zhang   } else if (isstring) { /* not test it yet! */
454b9ad928SBarry Smith     ierr = PetscViewerStringSPrintf(viewer," Redundant solver preconditioner");CHKERRQ(ierr);
464b9ad928SBarry Smith     ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr);
474b9ad928SBarry Smith     if (!rank) {
484b9ad928SBarry Smith       ierr = PCView(red->pc,sviewer);CHKERRQ(ierr);
494b9ad928SBarry Smith     }
504b9ad928SBarry Smith     ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr);
514b9ad928SBarry Smith   } else {
5279a5c55eSBarry Smith     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PC redundant",((PetscObject)viewer)->type_name);
534b9ad928SBarry Smith   }
544b9ad928SBarry Smith   PetscFunctionReturn(0);
554b9ad928SBarry Smith }
564b9ad928SBarry Smith 
57b9147fbbSdalcinl #include "include/private/matimpl.h"        /*I "petscmat.h" I*/
584b9ad928SBarry Smith #undef __FUNCT__
594b9ad928SBarry Smith #define __FUNCT__ "PCSetUp_Redundant"
606849ba73SBarry Smith static PetscErrorCode PCSetUp_Redundant(PC pc)
614b9ad928SBarry Smith {
624b9ad928SBarry Smith   PC_Redundant   *red = (PC_Redundant*)pc->data;
63dfbe8321SBarry Smith   PetscErrorCode ierr;
6445fc02eaSBarry Smith   PetscInt       mstart,mend,mlocal,m,mlocal_sub,rstart_sub,rend_sub,mloc_sub;
6513f74950SBarry Smith   PetscMPIInt    size;
664b9ad928SBarry Smith   MatReuse       reuse = MAT_INITIAL_MATRIX;
674b9ad928SBarry Smith   MatStructure   str = DIFFERENT_NONZERO_PATTERN;
680d7810c8SBarry Smith   MPI_Comm       comm = pc->comm,subcomm;
6923ce1328SBarry Smith   Vec            vec;
703f457be1SHong Zhang   PetscMPIInt    subsize,subrank;
711fbd8f88SHong Zhang   const char     *prefix;
723f457be1SHong Zhang 
734b9ad928SBarry Smith   PetscFunctionBegin;
7423ce1328SBarry Smith   ierr = MatGetVecs(pc->pmat,&vec,0);CHKERRQ(ierr);
7523ce1328SBarry Smith   ierr = VecGetSize(vec,&m);CHKERRQ(ierr);
761fbd8f88SHong Zhang 
774b9ad928SBarry Smith   if (!pc->setupcalled) {
781fbd8f88SHong Zhang     ierr = PetscSubcommCreate(comm,red->nsubcomm,&red->psubcomm);CHKERRQ(ierr);
791fbd8f88SHong Zhang     ierr = PetscLogObjectMemory(pc,sizeof(PetscSubcomm));CHKERRQ(ierr);
801fbd8f88SHong Zhang 
811fbd8f88SHong Zhang     /* create a new PC that processors in each subcomm have copy of */
820d7810c8SBarry Smith     subcomm = red->psubcomm->comm;
831fbd8f88SHong Zhang     ierr = PCCreate(subcomm,&red->pc);CHKERRQ(ierr);
841fbd8f88SHong Zhang     ierr = PCSetType(red->pc,PCLU);CHKERRQ(ierr);
851fbd8f88SHong Zhang     ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
861fbd8f88SHong Zhang     ierr = PCSetOptionsPrefix(red->pc,prefix);CHKERRQ(ierr);
871fbd8f88SHong Zhang     ierr = PCAppendOptionsPrefix(red->pc,"redundant_");CHKERRQ(ierr);
881fbd8f88SHong Zhang 
893f457be1SHong Zhang     /* create working vectors xsub/ysub and xdup/ydup */
9023ce1328SBarry Smith     ierr = VecGetLocalSize(vec,&mlocal);CHKERRQ(ierr);
913f457be1SHong Zhang     ierr = VecGetOwnershipRange(vec,&mstart,&mend);CHKERRQ(ierr);
924b9ad928SBarry Smith 
933f457be1SHong Zhang     /* get local size of xsub/ysub */
941fbd8f88SHong Zhang     ierr = MPI_Comm_size(subcomm,&subsize);CHKERRQ(ierr);
951fbd8f88SHong Zhang     ierr = MPI_Comm_rank(subcomm,&subrank);CHKERRQ(ierr);
961fbd8f88SHong Zhang     rstart_sub = pc->pmat->rmap.range[red->psubcomm->n*subrank]; /* rstart in xsub/ysub */
973f457be1SHong Zhang     if (subrank+1 < subsize){
981fbd8f88SHong Zhang       rend_sub = pc->pmat->rmap.range[red->psubcomm->n*(subrank+1)];
993f457be1SHong Zhang     } else {
1003f457be1SHong Zhang       rend_sub = m;
1013f457be1SHong Zhang     }
1023f457be1SHong Zhang     mloc_sub = rend_sub - rstart_sub;
1031fbd8f88SHong Zhang     ierr = VecCreateMPI(subcomm,mloc_sub,PETSC_DECIDE,&red->ysub);CHKERRQ(ierr);
1043f457be1SHong Zhang     /* create xsub with empty local arrays, because xdup's arrays will be placed into it */
1051fbd8f88SHong Zhang     ierr = VecCreateMPIWithArray(subcomm,mloc_sub,PETSC_DECIDE,PETSC_NULL,&red->xsub);CHKERRQ(ierr);
1063f457be1SHong Zhang 
1073f457be1SHong Zhang     /* create xdup and ydup. ydup has empty local arrays because ysub's arrays will be place into it.
1083f457be1SHong Zhang        Note: we use communicator dupcomm, not pc->comm! */
1091fbd8f88SHong Zhang     ierr = VecCreateMPI(red->psubcomm->dupparent,mloc_sub,PETSC_DECIDE,&red->xdup);CHKERRQ(ierr);
1101fbd8f88SHong Zhang     ierr = VecCreateMPIWithArray(red->psubcomm->dupparent,mloc_sub,PETSC_DECIDE,PETSC_NULL,&red->ydup);CHKERRQ(ierr);
1113f457be1SHong Zhang 
1123f457be1SHong Zhang     /* create vec scatters */
1133f457be1SHong Zhang     if (!red->scatterin){
1143f457be1SHong Zhang       IS       is1,is2;
1153f457be1SHong Zhang       PetscInt *idx1,*idx2,i,j,k;
11645fc02eaSBarry Smith 
1171fbd8f88SHong Zhang       ierr = PetscMalloc(2*red->psubcomm->n*mlocal*sizeof(PetscInt),&idx1);CHKERRQ(ierr);
1181fbd8f88SHong Zhang       idx2 = idx1 + red->psubcomm->n*mlocal;
1193f457be1SHong Zhang       j = 0;
1201fbd8f88SHong Zhang       for (k=0; k<red->psubcomm->n; k++){
1213f457be1SHong Zhang         for (i=mstart; i<mend; i++){
1223f457be1SHong Zhang           idx1[j]   = i;
1233f457be1SHong Zhang           idx2[j++] = i + m*k;
1243f457be1SHong Zhang         }
1253f457be1SHong Zhang       }
1261fbd8f88SHong Zhang       ierr = ISCreateGeneral(comm,red->psubcomm->n*mlocal,idx1,&is1);CHKERRQ(ierr);
1271fbd8f88SHong Zhang       ierr = ISCreateGeneral(comm,red->psubcomm->n*mlocal,idx2,&is2);CHKERRQ(ierr);
1283f457be1SHong Zhang       ierr = VecScatterCreate(vec,is1,red->xdup,is2,&red->scatterin);CHKERRQ(ierr);
1293f457be1SHong Zhang       ierr = ISDestroy(is1);CHKERRQ(ierr);
1303f457be1SHong Zhang       ierr = ISDestroy(is2);CHKERRQ(ierr);
1313f457be1SHong Zhang 
1321fbd8f88SHong Zhang       ierr = ISCreateStride(comm,mlocal,mstart+ red->psubcomm->color*m,1,&is1);CHKERRQ(ierr);
1333f457be1SHong Zhang       ierr = ISCreateStride(comm,mlocal,mstart,1,&is2);CHKERRQ(ierr);
1343f457be1SHong Zhang       ierr = VecScatterCreate(red->xdup,is1,vec,is2,&red->scatterout);CHKERRQ(ierr);
1353f457be1SHong Zhang       ierr = ISDestroy(is1);CHKERRQ(ierr);
1363f457be1SHong Zhang       ierr = ISDestroy(is2);CHKERRQ(ierr);
1373f457be1SHong Zhang       ierr = PetscFree(idx1);CHKERRQ(ierr);
1384b9ad928SBarry Smith     }
1394b9ad928SBarry Smith   }
14023ce1328SBarry Smith   ierr = VecDestroy(vec);CHKERRQ(ierr);
1414b9ad928SBarry Smith 
1424b9ad928SBarry Smith   /* if pmatrix set by user is sequential then we do not need to gather the parallel matrix */
1433f457be1SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1444b9ad928SBarry Smith   if (size == 1) {
1454b9ad928SBarry Smith     red->useparallelmat = PETSC_FALSE;
1464b9ad928SBarry Smith   }
1474b9ad928SBarry Smith 
1484b9ad928SBarry Smith   if (red->useparallelmat) {
1494b9ad928SBarry Smith     if (pc->setupcalled == 1 && pc->flag == DIFFERENT_NONZERO_PATTERN) {
1504b9ad928SBarry Smith       /* destroy old matrices */
1514b9ad928SBarry Smith       if (red->pmats) {
152b3804887SHong Zhang         ierr = MatDestroy(red->pmats);CHKERRQ(ierr);
1534b9ad928SBarry Smith       }
1544b9ad928SBarry Smith     } else if (pc->setupcalled == 1) {
1554b9ad928SBarry Smith       reuse = MAT_REUSE_MATRIX;
1564b9ad928SBarry Smith       str   = SAME_NONZERO_PATTERN;
1574b9ad928SBarry Smith     }
1584b9ad928SBarry Smith 
1593f457be1SHong Zhang     /* grab the parallel matrix and put it into processors of a subcomminicator */
160f664ae05SHong Zhang     /*--------------------------------------------------------------------------*/
161f664ae05SHong Zhang     ierr = VecGetLocalSize(red->ysub,&mlocal_sub);CHKERRQ(ierr);
16269db28dcSHong Zhang     ierr = MatGetRedundantMatrix(pc->pmat,red->psubcomm->n,red->psubcomm->comm,mlocal_sub,reuse,&red->pmats);CHKERRQ(ierr);
16369db28dcSHong Zhang 
1643f457be1SHong Zhang     /* tell PC of the subcommunicator its operators */
165b3804887SHong Zhang     ierr = PCSetOperators(red->pc,red->pmats,red->pmats,str);CHKERRQ(ierr);
1664b9ad928SBarry Smith   } else {
1674b9ad928SBarry Smith     ierr = PCSetOperators(red->pc,pc->mat,pc->pmat,pc->flag);CHKERRQ(ierr);
1684b9ad928SBarry Smith   }
169*0c24e6a1SHong Zhang   if (pc->setfromoptionscalled){
1704b9ad928SBarry Smith     ierr = PCSetFromOptions(red->pc);CHKERRQ(ierr);
171*0c24e6a1SHong Zhang   }
1724b9ad928SBarry Smith   ierr = PCSetUp(red->pc);CHKERRQ(ierr);
1734b9ad928SBarry Smith   PetscFunctionReturn(0);
1744b9ad928SBarry Smith }
1754b9ad928SBarry Smith 
1764b9ad928SBarry Smith #undef __FUNCT__
1774b9ad928SBarry Smith #define __FUNCT__ "PCApply_Redundant"
1786849ba73SBarry Smith static PetscErrorCode PCApply_Redundant(PC pc,Vec x,Vec y)
1794b9ad928SBarry Smith {
1804b9ad928SBarry Smith   PC_Redundant   *red = (PC_Redundant*)pc->data;
181dfbe8321SBarry Smith   PetscErrorCode ierr;
1823f457be1SHong Zhang   PetscScalar    *array;
1834b9ad928SBarry Smith 
1844b9ad928SBarry Smith   PetscFunctionBegin;
1853f457be1SHong Zhang   /* scatter x to xdup */
1863f457be1SHong Zhang   ierr = VecScatterBegin(x,red->xdup,INSERT_VALUES,SCATTER_FORWARD,red->scatterin);CHKERRQ(ierr);
1873f457be1SHong Zhang   ierr = VecScatterEnd(x,red->xdup,INSERT_VALUES,SCATTER_FORWARD,red->scatterin);CHKERRQ(ierr);
1883f457be1SHong Zhang 
1893f457be1SHong Zhang   /* place xdup's local array into xsub */
1903f457be1SHong Zhang   ierr = VecGetArray(red->xdup,&array);CHKERRQ(ierr);
1913f457be1SHong Zhang   ierr = VecPlaceArray(red->xsub,(const PetscScalar*)array);CHKERRQ(ierr);
1924b9ad928SBarry Smith 
1934b9ad928SBarry Smith   /* apply preconditioner on each processor */
1943f457be1SHong Zhang   ierr = PCApply(red->pc,red->xsub,red->ysub);CHKERRQ(ierr);
1953f457be1SHong Zhang   ierr = VecResetArray(red->xsub);CHKERRQ(ierr);
1963f457be1SHong Zhang   ierr = VecRestoreArray(red->xdup,&array);CHKERRQ(ierr);
1974b9ad928SBarry Smith 
1983f457be1SHong Zhang   /* place ysub's local array into ydup */
1993f457be1SHong Zhang   ierr = VecGetArray(red->ysub,&array);CHKERRQ(ierr);
2003f457be1SHong Zhang   ierr = VecPlaceArray(red->ydup,(const PetscScalar*)array);CHKERRQ(ierr);
2013f457be1SHong Zhang 
2023f457be1SHong Zhang   /* scatter ydup to y */
2033f457be1SHong Zhang   ierr = VecScatterBegin(red->ydup,y,INSERT_VALUES,SCATTER_FORWARD,red->scatterout);CHKERRQ(ierr);
2043f457be1SHong Zhang   ierr = VecScatterEnd(red->ydup,y,INSERT_VALUES,SCATTER_FORWARD,red->scatterout);CHKERRQ(ierr);
2053f457be1SHong Zhang   ierr = VecResetArray(red->ydup);CHKERRQ(ierr);
2063f457be1SHong Zhang   ierr = VecRestoreArray(red->ysub,&array);CHKERRQ(ierr);
2074b9ad928SBarry Smith   PetscFunctionReturn(0);
2084b9ad928SBarry Smith }
2094b9ad928SBarry Smith 
2104b9ad928SBarry Smith #undef __FUNCT__
2114b9ad928SBarry Smith #define __FUNCT__ "PCDestroy_Redundant"
2126849ba73SBarry Smith static PetscErrorCode PCDestroy_Redundant(PC pc)
2134b9ad928SBarry Smith {
2144b9ad928SBarry Smith   PC_Redundant   *red = (PC_Redundant*)pc->data;
215dfbe8321SBarry Smith   PetscErrorCode ierr;
2164b9ad928SBarry Smith 
2174b9ad928SBarry Smith   PetscFunctionBegin;
2184b9ad928SBarry Smith   if (red->scatterin)  {ierr = VecScatterDestroy(red->scatterin);CHKERRQ(ierr);}
2194b9ad928SBarry Smith   if (red->scatterout) {ierr = VecScatterDestroy(red->scatterout);CHKERRQ(ierr);}
2203f457be1SHong Zhang   if (red->ysub)       {ierr = VecDestroy(red->ysub);CHKERRQ(ierr);}
2213f457be1SHong Zhang   if (red->xsub)       {ierr = VecDestroy(red->xsub);CHKERRQ(ierr);}
2223f457be1SHong Zhang   if (red->xdup)       {ierr = VecDestroy(red->xdup);CHKERRQ(ierr);}
2233f457be1SHong Zhang   if (red->ydup)       {ierr = VecDestroy(red->ydup);CHKERRQ(ierr);}
224b3804887SHong Zhang   if (red->pmats) {
225b3804887SHong Zhang     ierr = MatDestroy(red->pmats);CHKERRQ(ierr);
2263f457be1SHong Zhang   }
227e5a9bf91SBarry Smith   if (red->psubcomm) {ierr = PetscSubcommDestroy(red->psubcomm);CHKERRQ(ierr);}
228e5a9bf91SBarry Smith   if (red->pc) {ierr = PCDestroy(red->pc);CHKERRQ(ierr);}
2294b9ad928SBarry Smith   ierr = PetscFree(red);CHKERRQ(ierr);
2304b9ad928SBarry Smith   PetscFunctionReturn(0);
2314b9ad928SBarry Smith }
2324b9ad928SBarry Smith 
2334b9ad928SBarry Smith #undef __FUNCT__
2344b9ad928SBarry Smith #define __FUNCT__ "PCSetFromOptions_Redundant"
2356849ba73SBarry Smith static PetscErrorCode PCSetFromOptions_Redundant(PC pc)
2364b9ad928SBarry Smith {
237a98ce0f4SHong Zhang   PetscErrorCode ierr;
238a98ce0f4SHong Zhang   PC_Redundant   *red = (PC_Redundant*)pc->data;
239a98ce0f4SHong Zhang 
2404b9ad928SBarry Smith   PetscFunctionBegin;
241a98ce0f4SHong Zhang   ierr = PetscOptionsHead("Redundant options");CHKERRQ(ierr);
24209a6bc64SHong Zhang   ierr = PetscOptionsInt("-pc_redundant_number","Number of redundant pc","PCRedundantSetNumber",red->nsubcomm,&red->nsubcomm,0);CHKERRQ(ierr);
243a98ce0f4SHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
2444b9ad928SBarry Smith   PetscFunctionReturn(0);
2454b9ad928SBarry Smith }
2464b9ad928SBarry Smith 
2474b9ad928SBarry Smith EXTERN_C_BEGIN
2484b9ad928SBarry Smith #undef __FUNCT__
24909a6bc64SHong Zhang #define __FUNCT__ "PCRedundantSetNumber_Redundant"
25009a6bc64SHong Zhang PetscErrorCode PETSCKSP_DLLEXPORT PCRedundantSetNumber_Redundant(PC pc,PetscInt nreds)
25109a6bc64SHong Zhang {
25209a6bc64SHong Zhang   PC_Redundant *red = (PC_Redundant*)pc->data;
25309a6bc64SHong Zhang 
25409a6bc64SHong Zhang   PetscFunctionBegin;
25509a6bc64SHong Zhang   red->nsubcomm = nreds;
25609a6bc64SHong Zhang   PetscFunctionReturn(0);
25709a6bc64SHong Zhang }
25809a6bc64SHong Zhang EXTERN_C_END
25909a6bc64SHong Zhang 
26009a6bc64SHong Zhang #undef __FUNCT__
26109a6bc64SHong Zhang #define __FUNCT__ "PCRedundantSetNumber"
26209a6bc64SHong Zhang /*@
26309a6bc64SHong Zhang    PCRedundantSetNumber - Sets the number of redundant preconditioner contexts.
26409a6bc64SHong Zhang 
26509a6bc64SHong Zhang    Collective on PC
26609a6bc64SHong Zhang 
26709a6bc64SHong Zhang    Input Parameters:
26809a6bc64SHong Zhang +  pc - the preconditioner context
2699b21d695SBarry Smith -  nredundant - number of redundant preconditioner contexts; for example if you are using 64 MPI processes and
2709b21d695SBarry Smith                               use an nredundant of 4 there will be 4 parallel solves each on 16 = 64/4 processes.
27109a6bc64SHong Zhang 
27209a6bc64SHong Zhang    Level: advanced
27309a6bc64SHong Zhang 
27409a6bc64SHong Zhang .keywords: PC, redundant solve
27509a6bc64SHong Zhang @*/
27609a6bc64SHong Zhang PetscErrorCode PETSCKSP_DLLEXPORT PCRedundantSetNumber(PC pc,PetscInt nredundant)
27709a6bc64SHong Zhang {
27809a6bc64SHong Zhang   PetscErrorCode ierr,(*f)(PC,PetscInt);
27909a6bc64SHong Zhang 
28009a6bc64SHong Zhang   PetscFunctionBegin;
28109a6bc64SHong Zhang   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
28209a6bc64SHong Zhang   if (nredundant <= 0) SETERRQ1(PETSC_ERR_ARG_WRONG, "num of redundant pc %D must be positive",nredundant);
28309a6bc64SHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCRedundantSetNumber_C",(void (**)(void))&f);CHKERRQ(ierr);
28409a6bc64SHong Zhang   if (f) {
28509a6bc64SHong Zhang     ierr = (*f)(pc,nredundant);CHKERRQ(ierr);
28609a6bc64SHong Zhang   }
28709a6bc64SHong Zhang   PetscFunctionReturn(0);
28809a6bc64SHong Zhang }
28909a6bc64SHong Zhang 
29009a6bc64SHong Zhang EXTERN_C_BEGIN
29109a6bc64SHong Zhang #undef __FUNCT__
2924b9ad928SBarry Smith #define __FUNCT__ "PCRedundantSetScatter_Redundant"
293dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCRedundantSetScatter_Redundant(PC pc,VecScatter in,VecScatter out)
2944b9ad928SBarry Smith {
2954b9ad928SBarry Smith   PC_Redundant   *red = (PC_Redundant*)pc->data;
296dfbe8321SBarry Smith   PetscErrorCode ierr;
2974b9ad928SBarry Smith 
2984b9ad928SBarry Smith   PetscFunctionBegin;
2994b9ad928SBarry Smith   ierr = PetscObjectReference((PetscObject)in);CHKERRQ(ierr);
300c3122656SLisandro Dalcin   if (red->scatterin) { ierr = VecScatterDestroy(red->scatterin);CHKERRQ(ierr); }
301c3122656SLisandro Dalcin   red->scatterin  = in;
3024b9ad928SBarry Smith   ierr = PetscObjectReference((PetscObject)out);CHKERRQ(ierr);
303c3122656SLisandro Dalcin   if (red->scatterout) { ierr = VecScatterDestroy(red->scatterout);CHKERRQ(ierr); }
304c3122656SLisandro Dalcin   red->scatterout = out;
3054b9ad928SBarry Smith   PetscFunctionReturn(0);
3064b9ad928SBarry Smith }
3074b9ad928SBarry Smith EXTERN_C_END
3084b9ad928SBarry Smith 
3094b9ad928SBarry Smith #undef __FUNCT__
3104b9ad928SBarry Smith #define __FUNCT__ "PCRedundantSetScatter"
3114b9ad928SBarry Smith /*@
3124b9ad928SBarry Smith    PCRedundantSetScatter - Sets the scatter used to copy values into the
3134b9ad928SBarry Smith      redundant local solve and the scatter to move them back into the global
3144b9ad928SBarry Smith      vector.
3154b9ad928SBarry Smith 
3164b9ad928SBarry Smith    Collective on PC
3174b9ad928SBarry Smith 
3184b9ad928SBarry Smith    Input Parameters:
3194b9ad928SBarry Smith +  pc - the preconditioner context
3204b9ad928SBarry Smith .  in - the scatter to move the values in
3214b9ad928SBarry Smith -  out - the scatter to move them out
3224b9ad928SBarry Smith 
3234b9ad928SBarry Smith    Level: advanced
3244b9ad928SBarry Smith 
3254b9ad928SBarry Smith .keywords: PC, redundant solve
3264b9ad928SBarry Smith @*/
327dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCRedundantSetScatter(PC pc,VecScatter in,VecScatter out)
3284b9ad928SBarry Smith {
329dfbe8321SBarry Smith   PetscErrorCode ierr,(*f)(PC,VecScatter,VecScatter);
3304b9ad928SBarry Smith 
3314b9ad928SBarry Smith   PetscFunctionBegin;
3324482741eSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
3334482741eSBarry Smith   PetscValidHeaderSpecific(in,VEC_SCATTER_COOKIE,2);
3344482741eSBarry Smith   PetscValidHeaderSpecific(out,VEC_SCATTER_COOKIE,3);
3354b9ad928SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCRedundantSetScatter_C",(void (**)(void))&f);CHKERRQ(ierr);
3364b9ad928SBarry Smith   if (f) {
3374b9ad928SBarry Smith     ierr = (*f)(pc,in,out);CHKERRQ(ierr);
3384b9ad928SBarry Smith   }
3394b9ad928SBarry Smith   PetscFunctionReturn(0);
3404b9ad928SBarry Smith }
3414b9ad928SBarry Smith 
3424b9ad928SBarry Smith EXTERN_C_BEGIN
3434b9ad928SBarry Smith #undef __FUNCT__
3444b9ad928SBarry Smith #define __FUNCT__ "PCRedundantGetPC_Redundant"
345dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCRedundantGetPC_Redundant(PC pc,PC *innerpc)
3464b9ad928SBarry Smith {
3474b9ad928SBarry Smith   PC_Redundant *red = (PC_Redundant*)pc->data;
3484b9ad928SBarry Smith 
3494b9ad928SBarry Smith   PetscFunctionBegin;
3504b9ad928SBarry Smith   *innerpc = red->pc;
3514b9ad928SBarry Smith   PetscFunctionReturn(0);
3524b9ad928SBarry Smith }
3534b9ad928SBarry Smith EXTERN_C_END
3544b9ad928SBarry Smith 
3554b9ad928SBarry Smith #undef __FUNCT__
3564b9ad928SBarry Smith #define __FUNCT__ "PCRedundantGetPC"
3574b9ad928SBarry Smith /*@
3584b9ad928SBarry Smith    PCRedundantGetPC - Gets the sequential PC created by the redundant PC.
3594b9ad928SBarry Smith 
3604b9ad928SBarry Smith    Not Collective
3614b9ad928SBarry Smith 
3624b9ad928SBarry Smith    Input Parameter:
3634b9ad928SBarry Smith .  pc - the preconditioner context
3644b9ad928SBarry Smith 
3654b9ad928SBarry Smith    Output Parameter:
3664b9ad928SBarry Smith .  innerpc - the sequential PC
3674b9ad928SBarry Smith 
3684b9ad928SBarry Smith    Level: advanced
3694b9ad928SBarry Smith 
3704b9ad928SBarry Smith .keywords: PC, redundant solve
3714b9ad928SBarry Smith @*/
372dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCRedundantGetPC(PC pc,PC *innerpc)
3734b9ad928SBarry Smith {
374dfbe8321SBarry Smith   PetscErrorCode ierr,(*f)(PC,PC*);
3754b9ad928SBarry Smith 
3764b9ad928SBarry Smith   PetscFunctionBegin;
3774482741eSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
3784482741eSBarry Smith   PetscValidPointer(innerpc,2);
3794b9ad928SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCRedundantGetPC_C",(void (**)(void))&f);CHKERRQ(ierr);
3804b9ad928SBarry Smith   if (f) {
3814b9ad928SBarry Smith     ierr = (*f)(pc,innerpc);CHKERRQ(ierr);
3824b9ad928SBarry Smith   }
3834b9ad928SBarry Smith   PetscFunctionReturn(0);
3844b9ad928SBarry Smith }
3854b9ad928SBarry Smith 
3864b9ad928SBarry Smith EXTERN_C_BEGIN
3874b9ad928SBarry Smith #undef __FUNCT__
3884b9ad928SBarry Smith #define __FUNCT__ "PCRedundantGetOperators_Redundant"
389dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCRedundantGetOperators_Redundant(PC pc,Mat *mat,Mat *pmat)
3904b9ad928SBarry Smith {
3914b9ad928SBarry Smith   PC_Redundant *red = (PC_Redundant*)pc->data;
3924b9ad928SBarry Smith 
3934b9ad928SBarry Smith   PetscFunctionBegin;
394b3804887SHong Zhang   if (mat)  *mat  = red->pmats;
395b3804887SHong Zhang   if (pmat) *pmat = red->pmats;
3964b9ad928SBarry Smith   PetscFunctionReturn(0);
3974b9ad928SBarry Smith }
3984b9ad928SBarry Smith EXTERN_C_END
3994b9ad928SBarry Smith 
4004b9ad928SBarry Smith #undef __FUNCT__
4014b9ad928SBarry Smith #define __FUNCT__ "PCRedundantGetOperators"
4024b9ad928SBarry Smith /*@
4034b9ad928SBarry Smith    PCRedundantGetOperators - gets the sequential matrix and preconditioner matrix
4044b9ad928SBarry Smith 
4054b9ad928SBarry Smith    Not Collective
4064b9ad928SBarry Smith 
4074b9ad928SBarry Smith    Input Parameter:
4084b9ad928SBarry Smith .  pc - the preconditioner context
4094b9ad928SBarry Smith 
4104b9ad928SBarry Smith    Output Parameters:
4114b9ad928SBarry Smith +  mat - the matrix
4124b9ad928SBarry Smith -  pmat - the (possibly different) preconditioner matrix
4134b9ad928SBarry Smith 
4144b9ad928SBarry Smith    Level: advanced
4154b9ad928SBarry Smith 
4164b9ad928SBarry Smith .keywords: PC, redundant solve
4174b9ad928SBarry Smith @*/
418dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCRedundantGetOperators(PC pc,Mat *mat,Mat *pmat)
4194b9ad928SBarry Smith {
420dfbe8321SBarry Smith   PetscErrorCode ierr,(*f)(PC,Mat*,Mat*);
4214b9ad928SBarry Smith 
4224b9ad928SBarry Smith   PetscFunctionBegin;
4234482741eSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
4244482741eSBarry Smith   if (mat)  PetscValidPointer(mat,2);
4254482741eSBarry Smith   if (pmat) PetscValidPointer(pmat,3);
4264b9ad928SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCRedundantGetOperators_C",(void (**)(void))&f);CHKERRQ(ierr);
4274b9ad928SBarry Smith   if (f) {
4284b9ad928SBarry Smith     ierr = (*f)(pc,mat,pmat);CHKERRQ(ierr);
4294b9ad928SBarry Smith   }
4304b9ad928SBarry Smith   PetscFunctionReturn(0);
4314b9ad928SBarry Smith }
4324b9ad928SBarry Smith 
4334b9ad928SBarry Smith /* -------------------------------------------------------------------------------------*/
43437a17b4dSBarry Smith /*MC
435a98ce0f4SHong Zhang      PCREDUNDANT - Runs a preconditioner for the entire problem on subgroups of processors
43637a17b4dSBarry Smith 
43737a17b4dSBarry Smith      Options for the redundant preconditioners can be set with -redundant_pc_xxx
43837a17b4dSBarry Smith 
43909391456SBarry Smith   Options Database:
4409b21d695SBarry Smith .  -pc_redundant_number <n> - number of redundant solves, for example if you are using 64 MPI processes and
4419b21d695SBarry Smith                               use an n of 4 there will be 4 parallel solves each on 16 = 64/4 processes.
44209391456SBarry Smith 
44337a17b4dSBarry Smith    Level: intermediate
44437a17b4dSBarry Smith 
44537a17b4dSBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PCRedundantSetScatter(),
4469b21d695SBarry Smith            PCRedundantGetPC(), PCRedundantGetOperators(), PCRedundantSetNumber()
44737a17b4dSBarry Smith M*/
44837a17b4dSBarry Smith 
4494b9ad928SBarry Smith EXTERN_C_BEGIN
4504b9ad928SBarry Smith #undef __FUNCT__
4514b9ad928SBarry Smith #define __FUNCT__ "PCCreate_Redundant"
452dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_Redundant(PC pc)
4534b9ad928SBarry Smith {
454dfbe8321SBarry Smith   PetscErrorCode ierr;
4554b9ad928SBarry Smith   PC_Redundant   *red;
45669db28dcSHong Zhang   PetscMPIInt    size;
4573f457be1SHong Zhang 
4584b9ad928SBarry Smith   PetscFunctionBegin;
4594b9ad928SBarry Smith   ierr = PetscNew(PC_Redundant,&red);CHKERRQ(ierr);
46052e6d16bSBarry Smith   ierr = PetscLogObjectMemory(pc,sizeof(PC_Redundant));CHKERRQ(ierr);
46169db28dcSHong Zhang   ierr = MPI_Comm_size(pc->comm,&size);CHKERRQ(ierr);
46269db28dcSHong Zhang   red->nsubcomm       = size;
4634b9ad928SBarry Smith   red->useparallelmat = PETSC_TRUE;
4641fbd8f88SHong Zhang   pc->data            = (void*)red;
4654b9ad928SBarry Smith 
4664b9ad928SBarry Smith   pc->ops->apply           = PCApply_Redundant;
4674b9ad928SBarry Smith   pc->ops->applytranspose  = 0;
4684b9ad928SBarry Smith   pc->ops->setup           = PCSetUp_Redundant;
4694b9ad928SBarry Smith   pc->ops->destroy         = PCDestroy_Redundant;
4704b9ad928SBarry Smith   pc->ops->setfromoptions  = PCSetFromOptions_Redundant;
4714b9ad928SBarry Smith   pc->ops->view            = PCView_Redundant;
4724b9ad928SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCRedundantSetScatter_C","PCRedundantSetScatter_Redundant",
4734b9ad928SBarry Smith                     PCRedundantSetScatter_Redundant);CHKERRQ(ierr);
47409a6bc64SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCRedundantSetNumber_C","PCRedundantSetNumber_Redundant",
47509a6bc64SHong Zhang                     PCRedundantSetNumber_Redundant);CHKERRQ(ierr);
4764b9ad928SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCRedundantGetPC_C","PCRedundantGetPC_Redundant",
4774b9ad928SBarry Smith                     PCRedundantGetPC_Redundant);CHKERRQ(ierr);
4784b9ad928SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCRedundantGetOperators_C","PCRedundantGetOperators_Redundant",
4794b9ad928SBarry Smith                     PCRedundantGetOperators_Redundant);CHKERRQ(ierr);
4804b9ad928SBarry Smith   PetscFunctionReturn(0);
4814b9ad928SBarry Smith }
4824b9ad928SBarry Smith EXTERN_C_END
483