xref: /petsc/src/ksp/pc/impls/redundant/redundant.c (revision 03ccd0b4bdafe6d76058ae38ab61bbb17e413f32)
1dba47a55SKris Buschelman #define PETSCKSP_DLL
2dba47a55SKris Buschelman 
34b9ad928SBarry Smith /*
43f457be1SHong Zhang   This file defines a "solve the problem redundantly on each subgroup of processor" preconditioner.
54b9ad928SBarry Smith */
66356e834SBarry Smith #include "private/pcimpl.h"     /*I "petscpc.h" I*/
74b9ad928SBarry Smith #include "petscksp.h"
84b9ad928SBarry Smith 
94b9ad928SBarry Smith typedef struct {
103e065800SHong Zhang   KSP          ksp;
114b9ad928SBarry Smith   PC           pc;                   /* actual preconditioner used on each processor */
127adad957SLisandro Dalcin   Vec          xsub,ysub;            /* vectors of a subcommunicator to hold parallel vectors of ((PetscObject)pc)->comm */
133f457be1SHong Zhang   Vec          xdup,ydup;            /* parallel vector that congregates xsub or ysub facilitating vector scattering */
14b3804887SHong Zhang   Mat          pmats;                /* matrix and optional preconditioner matrix belong to a subcommunicator */
153f457be1SHong Zhang   VecScatter   scatterin,scatterout; /* scatter used to move all values to each processor group (subcommunicator) */
164b9ad928SBarry Smith   PetscTruth   useparallelmat;
17c540e29cSHong Zhang   PetscSubcomm psubcomm;
181fbd8f88SHong Zhang   PetscInt     nsubcomm;           /* num of data structure PetscSubcomm */
194b9ad928SBarry Smith } PC_Redundant;
204b9ad928SBarry Smith 
214b9ad928SBarry Smith #undef __FUNCT__
224b9ad928SBarry Smith #define __FUNCT__ "PCView_Redundant"
236849ba73SBarry Smith static PetscErrorCode PCView_Redundant(PC pc,PetscViewer viewer)
244b9ad928SBarry Smith {
254b9ad928SBarry Smith   PC_Redundant   *red = (PC_Redundant*)pc->data;
26dfbe8321SBarry Smith   PetscErrorCode ierr;
2713f74950SBarry Smith   PetscMPIInt    rank;
2832077d6dSBarry Smith   PetscTruth     iascii,isstring;
29*03ccd0b4SBarry Smith   PetscViewer    subviewer;
304b9ad928SBarry Smith 
314b9ad928SBarry Smith   PetscFunctionBegin;
327adad957SLisandro Dalcin   ierr = MPI_Comm_rank(((PetscObject)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*03ccd0b4SBarry Smith     if (!red->psubcomm) {
37*03ccd0b4SBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  Redundant preconditioner: Not yet setup\n");CHKERRQ(ierr);
38*03ccd0b4SBarry Smith     } else {
393e065800SHong Zhang       ierr = PetscViewerASCIIPrintf(viewer,"  Redundant preconditioner: First (color=0) of %D PCs follows\n",red->nsubcomm);CHKERRQ(ierr);
407adad957SLisandro Dalcin       ierr = PetscViewerGetSubcomm(viewer,((PetscObject)red->pc)->comm,&subviewer);CHKERRQ(ierr);
41*03ccd0b4SBarry Smith       if (red->psubcomm->color) { /* only view first redundant pc */
424b9ad928SBarry Smith 	ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
433e065800SHong Zhang 	ierr = KSPView(red->ksp,subviewer);CHKERRQ(ierr);
444b9ad928SBarry Smith 	ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
454b9ad928SBarry Smith       }
467adad957SLisandro Dalcin       ierr = PetscViewerRestoreSubcomm(viewer,((PetscObject)red->pc)->comm,&subviewer);CHKERRQ(ierr);
474b9ad928SBarry Smith     }
48*03ccd0b4SBarry Smith   } else if (isstring) {
49*03ccd0b4SBarry Smith     ierr = PetscViewerStringSPrintf(viewer," Redundant solver preconditioner");CHKERRQ(ierr);
504b9ad928SBarry Smith   } else {
5179a5c55eSBarry Smith     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for PC redundant",((PetscObject)viewer)->type_name);
524b9ad928SBarry Smith   }
534b9ad928SBarry Smith   PetscFunctionReturn(0);
544b9ad928SBarry Smith }
554b9ad928SBarry Smith 
567c4f633dSBarry Smith #include "private/matimpl.h"        /*I "petscmat.h" I*/
574b9ad928SBarry Smith #undef __FUNCT__
584b9ad928SBarry Smith #define __FUNCT__ "PCSetUp_Redundant"
596849ba73SBarry Smith static PetscErrorCode PCSetUp_Redundant(PC pc)
604b9ad928SBarry Smith {
614b9ad928SBarry Smith   PC_Redundant   *red = (PC_Redundant*)pc->data;
62dfbe8321SBarry Smith   PetscErrorCode ierr;
6345fc02eaSBarry Smith   PetscInt       mstart,mend,mlocal,m,mlocal_sub,rstart_sub,rend_sub,mloc_sub;
6413f74950SBarry Smith   PetscMPIInt    size;
654b9ad928SBarry Smith   MatReuse       reuse = MAT_INITIAL_MATRIX;
664b9ad928SBarry Smith   MatStructure   str = DIFFERENT_NONZERO_PATTERN;
677adad957SLisandro Dalcin   MPI_Comm       comm = ((PetscObject)pc)->comm,subcomm;
6823ce1328SBarry Smith   Vec            vec;
693f457be1SHong Zhang   PetscMPIInt    subsize,subrank;
701fbd8f88SHong Zhang   const char     *prefix;
713f457be1SHong Zhang 
724b9ad928SBarry Smith   PetscFunctionBegin;
7323ce1328SBarry Smith   ierr = MatGetVecs(pc->pmat,&vec,0);CHKERRQ(ierr);
7423ce1328SBarry Smith   ierr = VecGetSize(vec,&m);CHKERRQ(ierr);
751fbd8f88SHong Zhang 
764b9ad928SBarry Smith   if (!pc->setupcalled) {
775f06b7aaSBarry Smith     if (!red->psubcomm) {
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;
835f06b7aaSBarry Smith       ierr = KSPCreate(subcomm,&red->ksp);CHKERRQ(ierr);
845f06b7aaSBarry Smith       ierr = PetscObjectIncrementTabLevel((PetscObject)red->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
855f06b7aaSBarry Smith       ierr = PetscLogObjectParent(pc,red->ksp);CHKERRQ(ierr);
865f06b7aaSBarry Smith       ierr = KSPSetType(red->ksp,KSPPREONLY);CHKERRQ(ierr);
875f06b7aaSBarry Smith       ierr = KSPGetPC(red->ksp,&red->pc);CHKERRQ(ierr);
88cf52b8b1SHong Zhang       ierr = PCSetType(red->pc,PCLU);CHKERRQ(ierr);
89cf52b8b1SHong Zhang 
901fbd8f88SHong Zhang       ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
915f06b7aaSBarry Smith       ierr = KSPSetOptionsPrefix(red->ksp,prefix);CHKERRQ(ierr);
925f06b7aaSBarry Smith       ierr = KSPAppendOptionsPrefix(red->ksp,"redundant_");CHKERRQ(ierr);
935f06b7aaSBarry Smith     } else {
945f06b7aaSBarry Smith        subcomm = red->psubcomm->comm;
955f06b7aaSBarry Smith     }
961fbd8f88SHong Zhang 
973f457be1SHong Zhang     /* create working vectors xsub/ysub and xdup/ydup */
9823ce1328SBarry Smith     ierr = VecGetLocalSize(vec,&mlocal);CHKERRQ(ierr);
993f457be1SHong Zhang     ierr = VecGetOwnershipRange(vec,&mstart,&mend);CHKERRQ(ierr);
1004b9ad928SBarry Smith 
1013f457be1SHong Zhang     /* get local size of xsub/ysub */
1021fbd8f88SHong Zhang     ierr = MPI_Comm_size(subcomm,&subsize);CHKERRQ(ierr);
1031fbd8f88SHong Zhang     ierr = MPI_Comm_rank(subcomm,&subrank);CHKERRQ(ierr);
104d0f46423SBarry Smith     rstart_sub = pc->pmat->rmap->range[red->psubcomm->n*subrank]; /* rstart in xsub/ysub */
1053f457be1SHong Zhang     if (subrank+1 < subsize){
106d0f46423SBarry Smith       rend_sub = pc->pmat->rmap->range[red->psubcomm->n*(subrank+1)];
1073f457be1SHong Zhang     } else {
1083f457be1SHong Zhang       rend_sub = m;
1093f457be1SHong Zhang     }
1103f457be1SHong Zhang     mloc_sub = rend_sub - rstart_sub;
1111fbd8f88SHong Zhang     ierr = VecCreateMPI(subcomm,mloc_sub,PETSC_DECIDE,&red->ysub);CHKERRQ(ierr);
1123f457be1SHong Zhang     /* create xsub with empty local arrays, because xdup's arrays will be placed into it */
1131fbd8f88SHong Zhang     ierr = VecCreateMPIWithArray(subcomm,mloc_sub,PETSC_DECIDE,PETSC_NULL,&red->xsub);CHKERRQ(ierr);
1143f457be1SHong Zhang 
1153f457be1SHong Zhang     /* create xdup and ydup. ydup has empty local arrays because ysub's arrays will be place into it.
1167adad957SLisandro Dalcin        Note: we use communicator dupcomm, not ((PetscObject)pc)->comm! */
1171fbd8f88SHong Zhang     ierr = VecCreateMPI(red->psubcomm->dupparent,mloc_sub,PETSC_DECIDE,&red->xdup);CHKERRQ(ierr);
1181fbd8f88SHong Zhang     ierr = VecCreateMPIWithArray(red->psubcomm->dupparent,mloc_sub,PETSC_DECIDE,PETSC_NULL,&red->ydup);CHKERRQ(ierr);
1193f457be1SHong Zhang 
1203f457be1SHong Zhang     /* create vec scatters */
1213f457be1SHong Zhang     if (!red->scatterin){
1223f457be1SHong Zhang       IS       is1,is2;
1233f457be1SHong Zhang       PetscInt *idx1,*idx2,i,j,k;
12445fc02eaSBarry Smith 
1251d79065fSBarry Smith       ierr = PetscMalloc2(red->psubcomm->n*mlocal,PetscInt,&idx1,red->psubcomm->n*mlocal,PetscInt,&idx2);CHKERRQ(ierr);
1263f457be1SHong Zhang       j = 0;
1271fbd8f88SHong Zhang       for (k=0; k<red->psubcomm->n; k++){
1283f457be1SHong Zhang         for (i=mstart; i<mend; i++){
1293f457be1SHong Zhang           idx1[j]   = i;
1303f457be1SHong Zhang           idx2[j++] = i + m*k;
1313f457be1SHong Zhang         }
1323f457be1SHong Zhang       }
1331fbd8f88SHong Zhang       ierr = ISCreateGeneral(comm,red->psubcomm->n*mlocal,idx1,&is1);CHKERRQ(ierr);
1341fbd8f88SHong Zhang       ierr = ISCreateGeneral(comm,red->psubcomm->n*mlocal,idx2,&is2);CHKERRQ(ierr);
1353f457be1SHong Zhang       ierr = VecScatterCreate(vec,is1,red->xdup,is2,&red->scatterin);CHKERRQ(ierr);
1363f457be1SHong Zhang       ierr = ISDestroy(is1);CHKERRQ(ierr);
1373f457be1SHong Zhang       ierr = ISDestroy(is2);CHKERRQ(ierr);
1383f457be1SHong Zhang 
1391fbd8f88SHong Zhang       ierr = ISCreateStride(comm,mlocal,mstart+ red->psubcomm->color*m,1,&is1);CHKERRQ(ierr);
1403f457be1SHong Zhang       ierr = ISCreateStride(comm,mlocal,mstart,1,&is2);CHKERRQ(ierr);
1413f457be1SHong Zhang       ierr = VecScatterCreate(red->xdup,is1,vec,is2,&red->scatterout);CHKERRQ(ierr);
1423f457be1SHong Zhang       ierr = ISDestroy(is1);CHKERRQ(ierr);
1433f457be1SHong Zhang       ierr = ISDestroy(is2);CHKERRQ(ierr);
1441d79065fSBarry Smith       ierr = PetscFree2(idx1,idx2);CHKERRQ(ierr);
1454b9ad928SBarry Smith     }
1464b9ad928SBarry Smith   }
14723ce1328SBarry Smith   ierr = VecDestroy(vec);CHKERRQ(ierr);
1484b9ad928SBarry Smith 
1494b9ad928SBarry Smith   /* if pmatrix set by user is sequential then we do not need to gather the parallel matrix */
1503f457be1SHong Zhang   ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1514b9ad928SBarry Smith   if (size == 1) {
1524b9ad928SBarry Smith     red->useparallelmat = PETSC_FALSE;
1534b9ad928SBarry Smith   }
1544b9ad928SBarry Smith 
1554b9ad928SBarry Smith   if (red->useparallelmat) {
1564b9ad928SBarry Smith     if (pc->setupcalled == 1 && pc->flag == DIFFERENT_NONZERO_PATTERN) {
1574b9ad928SBarry Smith       /* destroy old matrices */
1584b9ad928SBarry Smith       if (red->pmats) {
159b3804887SHong Zhang         ierr = MatDestroy(red->pmats);CHKERRQ(ierr);
1604b9ad928SBarry Smith       }
1614b9ad928SBarry Smith     } else if (pc->setupcalled == 1) {
1624b9ad928SBarry Smith       reuse = MAT_REUSE_MATRIX;
1634b9ad928SBarry Smith       str   = SAME_NONZERO_PATTERN;
1644b9ad928SBarry Smith     }
1654b9ad928SBarry Smith 
1663f457be1SHong Zhang     /* grab the parallel matrix and put it into processors of a subcomminicator */
167f664ae05SHong Zhang     /*--------------------------------------------------------------------------*/
168f664ae05SHong Zhang     ierr = VecGetLocalSize(red->ysub,&mlocal_sub);CHKERRQ(ierr);
16969db28dcSHong Zhang     ierr = MatGetRedundantMatrix(pc->pmat,red->psubcomm->n,red->psubcomm->comm,mlocal_sub,reuse,&red->pmats);CHKERRQ(ierr);
1703f457be1SHong Zhang     /* tell PC of the subcommunicator its operators */
17190f1c854SHong Zhang     ierr = KSPSetOperators(red->ksp,red->pmats,red->pmats,str);CHKERRQ(ierr);
1724b9ad928SBarry Smith   } else {
17390f1c854SHong Zhang     ierr = KSPSetOperators(red->ksp,pc->mat,pc->pmat,pc->flag);CHKERRQ(ierr);
1744b9ad928SBarry Smith   }
1750c24e6a1SHong Zhang   if (pc->setfromoptionscalled){
1763e065800SHong Zhang     ierr = KSPSetFromOptions(red->ksp);CHKERRQ(ierr);
1770c24e6a1SHong Zhang   }
1783e065800SHong Zhang   ierr = KSPSetUp(red->ksp);CHKERRQ(ierr);
1794b9ad928SBarry Smith   PetscFunctionReturn(0);
1804b9ad928SBarry Smith }
1814b9ad928SBarry Smith 
1824b9ad928SBarry Smith #undef __FUNCT__
1834b9ad928SBarry Smith #define __FUNCT__ "PCApply_Redundant"
1846849ba73SBarry Smith static PetscErrorCode PCApply_Redundant(PC pc,Vec x,Vec y)
1854b9ad928SBarry Smith {
1864b9ad928SBarry Smith   PC_Redundant   *red = (PC_Redundant*)pc->data;
187dfbe8321SBarry Smith   PetscErrorCode ierr;
1883f457be1SHong Zhang   PetscScalar    *array;
1894b9ad928SBarry Smith 
1904b9ad928SBarry Smith   PetscFunctionBegin;
1913f457be1SHong Zhang   /* scatter x to xdup */
192ca9f406cSSatish Balay   ierr = VecScatterBegin(red->scatterin,x,red->xdup,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
193ca9f406cSSatish Balay   ierr = VecScatterEnd(red->scatterin,x,red->xdup,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
1943f457be1SHong Zhang 
1953f457be1SHong Zhang   /* place xdup's local array into xsub */
1963f457be1SHong Zhang   ierr = VecGetArray(red->xdup,&array);CHKERRQ(ierr);
1973f457be1SHong Zhang   ierr = VecPlaceArray(red->xsub,(const PetscScalar*)array);CHKERRQ(ierr);
1984b9ad928SBarry Smith 
1994b9ad928SBarry Smith   /* apply preconditioner on each processor */
2003f457be1SHong Zhang   ierr = PCApply(red->pc,red->xsub,red->ysub);CHKERRQ(ierr);
2013f457be1SHong Zhang   ierr = VecResetArray(red->xsub);CHKERRQ(ierr);
2023f457be1SHong Zhang   ierr = VecRestoreArray(red->xdup,&array);CHKERRQ(ierr);
2034b9ad928SBarry Smith 
2043f457be1SHong Zhang   /* place ysub's local array into ydup */
2053f457be1SHong Zhang   ierr = VecGetArray(red->ysub,&array);CHKERRQ(ierr);
2063f457be1SHong Zhang   ierr = VecPlaceArray(red->ydup,(const PetscScalar*)array);CHKERRQ(ierr);
2073f457be1SHong Zhang 
2083f457be1SHong Zhang   /* scatter ydup to y */
209ca9f406cSSatish Balay   ierr = VecScatterBegin(red->scatterout,red->ydup,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
210ca9f406cSSatish Balay   ierr = VecScatterEnd(red->scatterout,red->ydup,y,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr);
2113f457be1SHong Zhang   ierr = VecResetArray(red->ydup);CHKERRQ(ierr);
2123f457be1SHong Zhang   ierr = VecRestoreArray(red->ysub,&array);CHKERRQ(ierr);
2134b9ad928SBarry Smith   PetscFunctionReturn(0);
2144b9ad928SBarry Smith }
2154b9ad928SBarry Smith 
2164b9ad928SBarry Smith #undef __FUNCT__
2174b9ad928SBarry Smith #define __FUNCT__ "PCDestroy_Redundant"
2186849ba73SBarry Smith static PetscErrorCode PCDestroy_Redundant(PC pc)
2194b9ad928SBarry Smith {
2204b9ad928SBarry Smith   PC_Redundant   *red = (PC_Redundant*)pc->data;
221dfbe8321SBarry Smith   PetscErrorCode ierr;
2224b9ad928SBarry Smith 
2234b9ad928SBarry Smith   PetscFunctionBegin;
2244b9ad928SBarry Smith   if (red->scatterin)  {ierr = VecScatterDestroy(red->scatterin);CHKERRQ(ierr);}
2254b9ad928SBarry Smith   if (red->scatterout) {ierr = VecScatterDestroy(red->scatterout);CHKERRQ(ierr);}
2263f457be1SHong Zhang   if (red->ysub)       {ierr = VecDestroy(red->ysub);CHKERRQ(ierr);}
2273f457be1SHong Zhang   if (red->xsub)       {ierr = VecDestroy(red->xsub);CHKERRQ(ierr);}
2283f457be1SHong Zhang   if (red->xdup)       {ierr = VecDestroy(red->xdup);CHKERRQ(ierr);}
2293f457be1SHong Zhang   if (red->ydup)       {ierr = VecDestroy(red->ydup);CHKERRQ(ierr);}
230b3804887SHong Zhang   if (red->pmats) {
231b3804887SHong Zhang     ierr = MatDestroy(red->pmats);CHKERRQ(ierr);
2323f457be1SHong Zhang   }
233e5a9bf91SBarry Smith   if (red->psubcomm) {ierr = PetscSubcommDestroy(red->psubcomm);CHKERRQ(ierr);}
2343e065800SHong Zhang   if (red->ksp) {ierr = KSPDestroy(red->ksp);CHKERRQ(ierr);}
2354b9ad928SBarry Smith   ierr = PetscFree(red);CHKERRQ(ierr);
2364b9ad928SBarry Smith   PetscFunctionReturn(0);
2374b9ad928SBarry Smith }
2384b9ad928SBarry Smith 
2394b9ad928SBarry Smith #undef __FUNCT__
2404b9ad928SBarry Smith #define __FUNCT__ "PCSetFromOptions_Redundant"
2416849ba73SBarry Smith static PetscErrorCode PCSetFromOptions_Redundant(PC pc)
2424b9ad928SBarry Smith {
243a98ce0f4SHong Zhang   PetscErrorCode ierr;
244a98ce0f4SHong Zhang   PC_Redundant   *red = (PC_Redundant*)pc->data;
245a98ce0f4SHong Zhang 
2464b9ad928SBarry Smith   PetscFunctionBegin;
247a98ce0f4SHong Zhang   ierr = PetscOptionsHead("Redundant options");CHKERRQ(ierr);
24809a6bc64SHong Zhang   ierr = PetscOptionsInt("-pc_redundant_number","Number of redundant pc","PCRedundantSetNumber",red->nsubcomm,&red->nsubcomm,0);CHKERRQ(ierr);
249a98ce0f4SHong Zhang   ierr = PetscOptionsTail();CHKERRQ(ierr);
2504b9ad928SBarry Smith   PetscFunctionReturn(0);
2514b9ad928SBarry Smith }
2524b9ad928SBarry Smith 
2534b9ad928SBarry Smith EXTERN_C_BEGIN
2544b9ad928SBarry Smith #undef __FUNCT__
25509a6bc64SHong Zhang #define __FUNCT__ "PCRedundantSetNumber_Redundant"
25609a6bc64SHong Zhang PetscErrorCode PETSCKSP_DLLEXPORT PCRedundantSetNumber_Redundant(PC pc,PetscInt nreds)
25709a6bc64SHong Zhang {
25809a6bc64SHong Zhang   PC_Redundant *red = (PC_Redundant*)pc->data;
25909a6bc64SHong Zhang 
26009a6bc64SHong Zhang   PetscFunctionBegin;
26109a6bc64SHong Zhang   red->nsubcomm = nreds;
26209a6bc64SHong Zhang   PetscFunctionReturn(0);
26309a6bc64SHong Zhang }
26409a6bc64SHong Zhang EXTERN_C_END
26509a6bc64SHong Zhang 
26609a6bc64SHong Zhang #undef __FUNCT__
26709a6bc64SHong Zhang #define __FUNCT__ "PCRedundantSetNumber"
26809a6bc64SHong Zhang /*@
26909a6bc64SHong Zhang    PCRedundantSetNumber - Sets the number of redundant preconditioner contexts.
27009a6bc64SHong Zhang 
27109a6bc64SHong Zhang    Collective on PC
27209a6bc64SHong Zhang 
27309a6bc64SHong Zhang    Input Parameters:
27409a6bc64SHong Zhang +  pc - the preconditioner context
2759b21d695SBarry Smith -  nredundant - number of redundant preconditioner contexts; for example if you are using 64 MPI processes and
2769b21d695SBarry Smith                               use an nredundant of 4 there will be 4 parallel solves each on 16 = 64/4 processes.
27709a6bc64SHong Zhang 
27809a6bc64SHong Zhang    Level: advanced
27909a6bc64SHong Zhang 
28009a6bc64SHong Zhang .keywords: PC, redundant solve
28109a6bc64SHong Zhang @*/
28209a6bc64SHong Zhang PetscErrorCode PETSCKSP_DLLEXPORT PCRedundantSetNumber(PC pc,PetscInt nredundant)
28309a6bc64SHong Zhang {
28409a6bc64SHong Zhang   PetscErrorCode ierr,(*f)(PC,PetscInt);
28509a6bc64SHong Zhang 
28609a6bc64SHong Zhang   PetscFunctionBegin;
28709a6bc64SHong Zhang   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
28809a6bc64SHong Zhang   if (nredundant <= 0) SETERRQ1(PETSC_ERR_ARG_WRONG, "num of redundant pc %D must be positive",nredundant);
28909a6bc64SHong Zhang   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCRedundantSetNumber_C",(void (**)(void))&f);CHKERRQ(ierr);
29009a6bc64SHong Zhang   if (f) {
29109a6bc64SHong Zhang     ierr = (*f)(pc,nredundant);CHKERRQ(ierr);
29209a6bc64SHong Zhang   }
29309a6bc64SHong Zhang   PetscFunctionReturn(0);
29409a6bc64SHong Zhang }
29509a6bc64SHong Zhang 
29609a6bc64SHong Zhang EXTERN_C_BEGIN
29709a6bc64SHong Zhang #undef __FUNCT__
2984b9ad928SBarry Smith #define __FUNCT__ "PCRedundantSetScatter_Redundant"
299dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCRedundantSetScatter_Redundant(PC pc,VecScatter in,VecScatter out)
3004b9ad928SBarry Smith {
3014b9ad928SBarry Smith   PC_Redundant   *red = (PC_Redundant*)pc->data;
302dfbe8321SBarry Smith   PetscErrorCode ierr;
3034b9ad928SBarry Smith 
3044b9ad928SBarry Smith   PetscFunctionBegin;
3054b9ad928SBarry Smith   ierr = PetscObjectReference((PetscObject)in);CHKERRQ(ierr);
306c3122656SLisandro Dalcin   if (red->scatterin) { ierr = VecScatterDestroy(red->scatterin);CHKERRQ(ierr); }
307c3122656SLisandro Dalcin   red->scatterin  = in;
3084b9ad928SBarry Smith   ierr = PetscObjectReference((PetscObject)out);CHKERRQ(ierr);
309c3122656SLisandro Dalcin   if (red->scatterout) { ierr = VecScatterDestroy(red->scatterout);CHKERRQ(ierr); }
310c3122656SLisandro Dalcin   red->scatterout = out;
3114b9ad928SBarry Smith   PetscFunctionReturn(0);
3124b9ad928SBarry Smith }
3134b9ad928SBarry Smith EXTERN_C_END
3144b9ad928SBarry Smith 
3154b9ad928SBarry Smith #undef __FUNCT__
3164b9ad928SBarry Smith #define __FUNCT__ "PCRedundantSetScatter"
3174b9ad928SBarry Smith /*@
3184b9ad928SBarry Smith    PCRedundantSetScatter - Sets the scatter used to copy values into the
3194b9ad928SBarry Smith      redundant local solve and the scatter to move them back into the global
3204b9ad928SBarry Smith      vector.
3214b9ad928SBarry Smith 
3224b9ad928SBarry Smith    Collective on PC
3234b9ad928SBarry Smith 
3244b9ad928SBarry Smith    Input Parameters:
3254b9ad928SBarry Smith +  pc - the preconditioner context
3264b9ad928SBarry Smith .  in - the scatter to move the values in
3274b9ad928SBarry Smith -  out - the scatter to move them out
3284b9ad928SBarry Smith 
3294b9ad928SBarry Smith    Level: advanced
3304b9ad928SBarry Smith 
3314b9ad928SBarry Smith .keywords: PC, redundant solve
3324b9ad928SBarry Smith @*/
333dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCRedundantSetScatter(PC pc,VecScatter in,VecScatter out)
3344b9ad928SBarry Smith {
335dfbe8321SBarry Smith   PetscErrorCode ierr,(*f)(PC,VecScatter,VecScatter);
3364b9ad928SBarry Smith 
3374b9ad928SBarry Smith   PetscFunctionBegin;
3384482741eSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
3394482741eSBarry Smith   PetscValidHeaderSpecific(in,VEC_SCATTER_COOKIE,2);
3404482741eSBarry Smith   PetscValidHeaderSpecific(out,VEC_SCATTER_COOKIE,3);
3414b9ad928SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCRedundantSetScatter_C",(void (**)(void))&f);CHKERRQ(ierr);
3424b9ad928SBarry Smith   if (f) {
3434b9ad928SBarry Smith     ierr = (*f)(pc,in,out);CHKERRQ(ierr);
3444b9ad928SBarry Smith   }
3454b9ad928SBarry Smith   PetscFunctionReturn(0);
3464b9ad928SBarry Smith }
3474b9ad928SBarry Smith 
3484b9ad928SBarry Smith EXTERN_C_BEGIN
3494b9ad928SBarry Smith #undef __FUNCT__
3504b9ad928SBarry Smith #define __FUNCT__ "PCRedundantGetPC_Redundant"
351dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCRedundantGetPC_Redundant(PC pc,PC *innerpc)
3524b9ad928SBarry Smith {
3535f06b7aaSBarry Smith   PetscErrorCode ierr;
3544b9ad928SBarry Smith   PC_Redundant   *red = (PC_Redundant*)pc->data;
3555f06b7aaSBarry Smith   MPI_Comm       comm,subcomm;
3565f06b7aaSBarry Smith   const char     *prefix;
3574b9ad928SBarry Smith 
3584b9ad928SBarry Smith   PetscFunctionBegin;
3595f06b7aaSBarry Smith   if (!red->psubcomm) {
3605f06b7aaSBarry Smith     ierr = PetscObjectGetComm((PetscObject)pc,&comm);CHKERRQ(ierr);
3615f06b7aaSBarry Smith     ierr = PetscSubcommCreate(comm,red->nsubcomm,&red->psubcomm);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;
3665f06b7aaSBarry Smith     ierr = KSPCreate(subcomm,&red->ksp);CHKERRQ(ierr);
3675f06b7aaSBarry Smith     ierr = PetscObjectIncrementTabLevel((PetscObject)red->ksp,(PetscObject)pc,1);CHKERRQ(ierr);
3685f06b7aaSBarry Smith     ierr = PetscLogObjectParent(pc,red->ksp);CHKERRQ(ierr);
3695f06b7aaSBarry Smith     ierr = KSPSetType(red->ksp,KSPPREONLY);CHKERRQ(ierr);
3705f06b7aaSBarry Smith     ierr = KSPGetPC(red->ksp,&red->pc);CHKERRQ(ierr);
3715f06b7aaSBarry Smith     ierr = PCSetType(red->pc,PCLU);CHKERRQ(ierr);
3725f06b7aaSBarry Smith 
3735f06b7aaSBarry Smith     ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr);
3745f06b7aaSBarry Smith     ierr = KSPSetOptionsPrefix(red->ksp,prefix);CHKERRQ(ierr);
3755f06b7aaSBarry Smith     ierr = KSPAppendOptionsPrefix(red->ksp,"redundant_");CHKERRQ(ierr);
3765f06b7aaSBarry Smith   }
3775f06b7aaSBarry Smith 
3785f06b7aaSBarry Smith   ierr = KSPGetPC(red->ksp,innerpc);CHKERRQ(ierr);
3794b9ad928SBarry Smith   PetscFunctionReturn(0);
3804b9ad928SBarry Smith }
3814b9ad928SBarry Smith EXTERN_C_END
3824b9ad928SBarry Smith 
3834b9ad928SBarry Smith #undef __FUNCT__
3844b9ad928SBarry Smith #define __FUNCT__ "PCRedundantGetPC"
3854b9ad928SBarry Smith /*@
3864b9ad928SBarry Smith    PCRedundantGetPC - Gets the sequential PC 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:
3944b9ad928SBarry Smith .  innerpc - the sequential PC
3954b9ad928SBarry Smith 
3964b9ad928SBarry Smith    Level: advanced
3974b9ad928SBarry Smith 
3984b9ad928SBarry Smith .keywords: PC, redundant solve
3994b9ad928SBarry Smith @*/
400dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCRedundantGetPC(PC pc,PC *innerpc)
4014b9ad928SBarry Smith {
402dfbe8321SBarry Smith   PetscErrorCode ierr,(*f)(PC,PC*);
4034b9ad928SBarry Smith 
4044b9ad928SBarry Smith   PetscFunctionBegin;
4054482741eSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
4064482741eSBarry Smith   PetscValidPointer(innerpc,2);
4074b9ad928SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCRedundantGetPC_C",(void (**)(void))&f);CHKERRQ(ierr);
4084b9ad928SBarry Smith   if (f) {
4094b9ad928SBarry Smith     ierr = (*f)(pc,innerpc);CHKERRQ(ierr);
4104b9ad928SBarry Smith   }
4114b9ad928SBarry Smith   PetscFunctionReturn(0);
4124b9ad928SBarry Smith }
4134b9ad928SBarry Smith 
4144b9ad928SBarry Smith EXTERN_C_BEGIN
4154b9ad928SBarry Smith #undef __FUNCT__
4164b9ad928SBarry Smith #define __FUNCT__ "PCRedundantGetOperators_Redundant"
417dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCRedundantGetOperators_Redundant(PC pc,Mat *mat,Mat *pmat)
4184b9ad928SBarry Smith {
4194b9ad928SBarry Smith   PC_Redundant *red = (PC_Redundant*)pc->data;
4204b9ad928SBarry Smith 
4214b9ad928SBarry Smith   PetscFunctionBegin;
422b3804887SHong Zhang   if (mat)  *mat  = red->pmats;
423b3804887SHong Zhang   if (pmat) *pmat = red->pmats;
4244b9ad928SBarry Smith   PetscFunctionReturn(0);
4254b9ad928SBarry Smith }
4264b9ad928SBarry Smith EXTERN_C_END
4274b9ad928SBarry Smith 
4284b9ad928SBarry Smith #undef __FUNCT__
4294b9ad928SBarry Smith #define __FUNCT__ "PCRedundantGetOperators"
4304b9ad928SBarry Smith /*@
4314b9ad928SBarry Smith    PCRedundantGetOperators - gets the sequential matrix and preconditioner matrix
4324b9ad928SBarry Smith 
4334b9ad928SBarry Smith    Not Collective
4344b9ad928SBarry Smith 
4354b9ad928SBarry Smith    Input Parameter:
4364b9ad928SBarry Smith .  pc - the preconditioner context
4374b9ad928SBarry Smith 
4384b9ad928SBarry Smith    Output Parameters:
4394b9ad928SBarry Smith +  mat - the matrix
4404b9ad928SBarry Smith -  pmat - the (possibly different) preconditioner matrix
4414b9ad928SBarry Smith 
4424b9ad928SBarry Smith    Level: advanced
4434b9ad928SBarry Smith 
4444b9ad928SBarry Smith .keywords: PC, redundant solve
4454b9ad928SBarry Smith @*/
446dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCRedundantGetOperators(PC pc,Mat *mat,Mat *pmat)
4474b9ad928SBarry Smith {
448dfbe8321SBarry Smith   PetscErrorCode ierr,(*f)(PC,Mat*,Mat*);
4494b9ad928SBarry Smith 
4504b9ad928SBarry Smith   PetscFunctionBegin;
4514482741eSBarry Smith   PetscValidHeaderSpecific(pc,PC_COOKIE,1);
4524482741eSBarry Smith   if (mat)  PetscValidPointer(mat,2);
4534482741eSBarry Smith   if (pmat) PetscValidPointer(pmat,3);
4544b9ad928SBarry Smith   ierr = PetscObjectQueryFunction((PetscObject)pc,"PCRedundantGetOperators_C",(void (**)(void))&f);CHKERRQ(ierr);
4554b9ad928SBarry Smith   if (f) {
4564b9ad928SBarry Smith     ierr = (*f)(pc,mat,pmat);CHKERRQ(ierr);
4574b9ad928SBarry Smith   }
4584b9ad928SBarry Smith   PetscFunctionReturn(0);
4594b9ad928SBarry Smith }
4604b9ad928SBarry Smith 
4614b9ad928SBarry Smith /* -------------------------------------------------------------------------------------*/
46237a17b4dSBarry Smith /*MC
463a98ce0f4SHong Zhang      PCREDUNDANT - Runs a preconditioner for the entire problem on subgroups of processors
46437a17b4dSBarry Smith 
46537a17b4dSBarry Smith      Options for the redundant preconditioners can be set with -redundant_pc_xxx
46637a17b4dSBarry Smith 
46709391456SBarry Smith   Options Database:
4689b21d695SBarry Smith .  -pc_redundant_number <n> - number of redundant solves, for example if you are using 64 MPI processes and
4699b21d695SBarry Smith                               use an n of 4 there will be 4 parallel solves each on 16 = 64/4 processes.
47009391456SBarry Smith 
47137a17b4dSBarry Smith    Level: intermediate
47237a17b4dSBarry Smith 
47337a17b4dSBarry Smith .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PCRedundantSetScatter(),
4749b21d695SBarry Smith            PCRedundantGetPC(), PCRedundantGetOperators(), PCRedundantSetNumber()
47537a17b4dSBarry Smith M*/
47637a17b4dSBarry Smith 
4774b9ad928SBarry Smith EXTERN_C_BEGIN
4784b9ad928SBarry Smith #undef __FUNCT__
4794b9ad928SBarry Smith #define __FUNCT__ "PCCreate_Redundant"
480dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_Redundant(PC pc)
4814b9ad928SBarry Smith {
482dfbe8321SBarry Smith   PetscErrorCode ierr;
4834b9ad928SBarry Smith   PC_Redundant   *red;
48469db28dcSHong Zhang   PetscMPIInt    size;
4853f457be1SHong Zhang 
4864b9ad928SBarry Smith   PetscFunctionBegin;
48738f2d2fdSLisandro Dalcin   ierr = PetscNewLog(pc,PC_Redundant,&red);CHKERRQ(ierr);
4887adad957SLisandro Dalcin   ierr = MPI_Comm_size(((PetscObject)pc)->comm,&size);CHKERRQ(ierr);
48969db28dcSHong Zhang   red->nsubcomm       = size;
4904b9ad928SBarry Smith   red->useparallelmat = PETSC_TRUE;
4911fbd8f88SHong Zhang   pc->data            = (void*)red;
4924b9ad928SBarry Smith 
4934b9ad928SBarry Smith   pc->ops->apply           = PCApply_Redundant;
4944b9ad928SBarry Smith   pc->ops->applytranspose  = 0;
4954b9ad928SBarry Smith   pc->ops->setup           = PCSetUp_Redundant;
4964b9ad928SBarry Smith   pc->ops->destroy         = PCDestroy_Redundant;
4974b9ad928SBarry Smith   pc->ops->setfromoptions  = PCSetFromOptions_Redundant;
4984b9ad928SBarry Smith   pc->ops->view            = PCView_Redundant;
4994b9ad928SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCRedundantSetScatter_C","PCRedundantSetScatter_Redundant",
5004b9ad928SBarry Smith                     PCRedundantSetScatter_Redundant);CHKERRQ(ierr);
50109a6bc64SHong Zhang   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCRedundantSetNumber_C","PCRedundantSetNumber_Redundant",
50209a6bc64SHong Zhang                     PCRedundantSetNumber_Redundant);CHKERRQ(ierr);
5034b9ad928SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCRedundantGetPC_C","PCRedundantGetPC_Redundant",
5044b9ad928SBarry Smith                     PCRedundantGetPC_Redundant);CHKERRQ(ierr);
5054b9ad928SBarry Smith   ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCRedundantGetOperators_C","PCRedundantGetOperators_Redundant",
5064b9ad928SBarry Smith                     PCRedundantGetOperators_Redundant);CHKERRQ(ierr);
5074b9ad928SBarry Smith   PetscFunctionReturn(0);
5084b9ad928SBarry Smith }
5094b9ad928SBarry Smith EXTERN_C_END
510