14b9ad928SBarry Smith /* 24b9ad928SBarry Smith This file defines a "solve the problem redundantly on each processor" preconditioner. 34b9ad928SBarry Smith 44b9ad928SBarry Smith */ 54b9ad928SBarry Smith #include "src/ksp/pc/pcimpl.h" /*I "petscpc.h" I*/ 64b9ad928SBarry Smith #include "petscksp.h" 74b9ad928SBarry Smith 84b9ad928SBarry Smith typedef struct { 94b9ad928SBarry Smith PC pc; /* actual preconditioner used on each processor */ 104b9ad928SBarry Smith Vec x,b; /* sequential vectors to hold parallel vectors */ 114b9ad928SBarry Smith Mat *pmats; /* matrix and optional preconditioner matrix */ 124b9ad928SBarry Smith VecScatter scatterin,scatterout; /* scatter used to move all values to each processor */ 134b9ad928SBarry Smith PetscTruth useparallelmat; 144b9ad928SBarry Smith } PC_Redundant; 154b9ad928SBarry Smith 164b9ad928SBarry Smith #undef __FUNCT__ 174b9ad928SBarry Smith #define __FUNCT__ "PCView_Redundant" 18*6849ba73SBarry Smith static PetscErrorCode PCView_Redundant(PC pc,PetscViewer viewer) 194b9ad928SBarry Smith { 204b9ad928SBarry Smith PC_Redundant *red = (PC_Redundant*)pc->data; 21dfbe8321SBarry Smith PetscErrorCode ierr; 22dfbe8321SBarry Smith int rank; 2332077d6dSBarry Smith PetscTruth iascii,isstring; 244b9ad928SBarry Smith PetscViewer sviewer; 254b9ad928SBarry Smith 264b9ad928SBarry Smith PetscFunctionBegin; 274b9ad928SBarry Smith ierr = MPI_Comm_rank(pc->comm,&rank);CHKERRQ(ierr); 2832077d6dSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 294b9ad928SBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr); 3032077d6dSBarry Smith if (iascii) { 314b9ad928SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Redundant solver preconditioner: Actual PC follows\n");CHKERRQ(ierr); 324b9ad928SBarry Smith ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr); 334b9ad928SBarry Smith if (!rank) { 344b9ad928SBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 354b9ad928SBarry Smith ierr = PCView(red->pc,sviewer);CHKERRQ(ierr); 364b9ad928SBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 374b9ad928SBarry Smith } 384b9ad928SBarry Smith ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr); 394b9ad928SBarry Smith } else if (isstring) { 404b9ad928SBarry Smith ierr = PetscViewerStringSPrintf(viewer," Redundant solver preconditioner");CHKERRQ(ierr); 414b9ad928SBarry Smith ierr = PetscViewerGetSingleton(viewer,&sviewer);CHKERRQ(ierr); 424b9ad928SBarry Smith if (!rank) { 434b9ad928SBarry Smith ierr = PCView(red->pc,sviewer);CHKERRQ(ierr); 444b9ad928SBarry Smith } 454b9ad928SBarry Smith ierr = PetscViewerRestoreSingleton(viewer,&sviewer);CHKERRQ(ierr); 464b9ad928SBarry Smith } else { 474b9ad928SBarry Smith SETERRQ1(1,"Viewer type %s not supported for PC redundant",((PetscObject)viewer)->type_name); 484b9ad928SBarry Smith } 494b9ad928SBarry Smith PetscFunctionReturn(0); 504b9ad928SBarry Smith } 514b9ad928SBarry Smith 524b9ad928SBarry Smith #undef __FUNCT__ 534b9ad928SBarry Smith #define __FUNCT__ "PCSetUp_Redundant" 54*6849ba73SBarry Smith static PetscErrorCode PCSetUp_Redundant(PC pc) 554b9ad928SBarry Smith { 564b9ad928SBarry Smith PC_Redundant *red = (PC_Redundant*)pc->data; 57dfbe8321SBarry Smith PetscErrorCode ierr; 58dfbe8321SBarry Smith int mstart,mlocal,m,size; 594b9ad928SBarry Smith IS isl; 604b9ad928SBarry Smith MatReuse reuse = MAT_INITIAL_MATRIX; 614b9ad928SBarry Smith MatStructure str = DIFFERENT_NONZERO_PATTERN; 624b9ad928SBarry Smith MPI_Comm comm; 6323ce1328SBarry Smith Vec vec; 644b9ad928SBarry Smith 654b9ad928SBarry Smith PetscFunctionBegin; 6623ce1328SBarry Smith ierr = MatGetVecs(pc->pmat,&vec,0);CHKERRQ(ierr); 674b9ad928SBarry Smith ierr = PCSetFromOptions(red->pc);CHKERRQ(ierr); 6823ce1328SBarry Smith ierr = VecGetSize(vec,&m);CHKERRQ(ierr); 694b9ad928SBarry Smith if (!pc->setupcalled) { 7023ce1328SBarry Smith ierr = VecGetLocalSize(vec,&mlocal);CHKERRQ(ierr); 714b9ad928SBarry Smith ierr = VecCreateSeq(PETSC_COMM_SELF,m,&red->x);CHKERRQ(ierr); 724b9ad928SBarry Smith ierr = VecDuplicate(red->x,&red->b);CHKERRQ(ierr); 734b9ad928SBarry Smith if (!red->scatterin) { 744b9ad928SBarry Smith 754b9ad928SBarry Smith /* 764b9ad928SBarry Smith Create the vectors and vector scatter to get the entire vector onto each processor 774b9ad928SBarry Smith */ 7823ce1328SBarry Smith ierr = VecGetOwnershipRange(vec,&mstart,PETSC_NULL);CHKERRQ(ierr); 7923ce1328SBarry Smith ierr = VecScatterCreate(vec,0,red->x,0,&red->scatterin);CHKERRQ(ierr); 804b9ad928SBarry Smith ierr = ISCreateStride(pc->comm,mlocal,mstart,1,&isl);CHKERRQ(ierr); 8123ce1328SBarry Smith ierr = VecScatterCreate(red->x,isl,vec,isl,&red->scatterout);CHKERRQ(ierr); 824b9ad928SBarry Smith ierr = ISDestroy(isl);CHKERRQ(ierr); 834b9ad928SBarry Smith } 844b9ad928SBarry Smith } 8523ce1328SBarry Smith ierr = VecDestroy(vec);CHKERRQ(ierr); 864b9ad928SBarry Smith 874b9ad928SBarry Smith /* if pmatrix set by user is sequential then we do not need to gather the parallel matrix*/ 884b9ad928SBarry Smith 894b9ad928SBarry Smith ierr = PetscObjectGetComm((PetscObject)pc->pmat,&comm);CHKERRQ(ierr); 904b9ad928SBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);CHKERRQ(ierr); 914b9ad928SBarry Smith if (size == 1) { 924b9ad928SBarry Smith red->useparallelmat = PETSC_FALSE; 934b9ad928SBarry Smith } 944b9ad928SBarry Smith 954b9ad928SBarry Smith if (red->useparallelmat) { 964b9ad928SBarry Smith if (pc->setupcalled == 1 && pc->flag == DIFFERENT_NONZERO_PATTERN) { 974b9ad928SBarry Smith /* destroy old matrices */ 984b9ad928SBarry Smith if (red->pmats) { 994b9ad928SBarry Smith ierr = MatDestroyMatrices(1,&red->pmats);CHKERRQ(ierr); 1004b9ad928SBarry Smith } 1014b9ad928SBarry Smith } else if (pc->setupcalled == 1) { 1024b9ad928SBarry Smith reuse = MAT_REUSE_MATRIX; 1034b9ad928SBarry Smith str = SAME_NONZERO_PATTERN; 1044b9ad928SBarry Smith } 1054b9ad928SBarry Smith 1064b9ad928SBarry Smith /* 1074b9ad928SBarry Smith grab the parallel matrix and put it on each processor 1084b9ad928SBarry Smith */ 1094b9ad928SBarry Smith ierr = ISCreateStride(PETSC_COMM_SELF,m,0,1,&isl);CHKERRQ(ierr); 1104b9ad928SBarry Smith ierr = MatGetSubMatrices(pc->pmat,1,&isl,&isl,reuse,&red->pmats);CHKERRQ(ierr); 1114b9ad928SBarry Smith ierr = ISDestroy(isl);CHKERRQ(ierr); 1124b9ad928SBarry Smith 1134b9ad928SBarry Smith /* tell sequential PC its operators */ 1144b9ad928SBarry Smith ierr = PCSetOperators(red->pc,red->pmats[0],red->pmats[0],str);CHKERRQ(ierr); 1154b9ad928SBarry Smith } else { 1164b9ad928SBarry Smith ierr = PCSetOperators(red->pc,pc->mat,pc->pmat,pc->flag);CHKERRQ(ierr); 1174b9ad928SBarry Smith } 1184b9ad928SBarry Smith ierr = PCSetFromOptions(red->pc);CHKERRQ(ierr); 1194b9ad928SBarry Smith ierr = PCSetUp(red->pc);CHKERRQ(ierr); 1204b9ad928SBarry Smith PetscFunctionReturn(0); 1214b9ad928SBarry Smith } 1224b9ad928SBarry Smith 1234b9ad928SBarry Smith 1244b9ad928SBarry Smith #undef __FUNCT__ 1254b9ad928SBarry Smith #define __FUNCT__ "PCApply_Redundant" 126*6849ba73SBarry Smith static PetscErrorCode PCApply_Redundant(PC pc,Vec x,Vec y) 1274b9ad928SBarry Smith { 1284b9ad928SBarry Smith PC_Redundant *red = (PC_Redundant*)pc->data; 129dfbe8321SBarry Smith PetscErrorCode ierr; 1304b9ad928SBarry Smith 1314b9ad928SBarry Smith PetscFunctionBegin; 1324b9ad928SBarry Smith /* move all values to each processor */ 1334b9ad928SBarry Smith ierr = VecScatterBegin(x,red->b,INSERT_VALUES,SCATTER_FORWARD,red->scatterin);CHKERRQ(ierr); 1344b9ad928SBarry Smith ierr = VecScatterEnd(x,red->b,INSERT_VALUES,SCATTER_FORWARD,red->scatterin);CHKERRQ(ierr); 1354b9ad928SBarry Smith 1364b9ad928SBarry Smith /* apply preconditioner on each processor */ 137d8fd42c4SBarry Smith ierr = PCApply(red->pc,red->b,red->x);CHKERRQ(ierr); 1384b9ad928SBarry Smith 1394b9ad928SBarry Smith /* move local part of values into y vector */ 1404b9ad928SBarry Smith ierr = VecScatterBegin(red->x,y,INSERT_VALUES,SCATTER_FORWARD,red->scatterout);CHKERRQ(ierr); 1414b9ad928SBarry Smith ierr = VecScatterEnd(red->x,y,INSERT_VALUES,SCATTER_FORWARD,red->scatterout);CHKERRQ(ierr); 1424b9ad928SBarry Smith PetscFunctionReturn(0); 1434b9ad928SBarry Smith } 1444b9ad928SBarry Smith 1454b9ad928SBarry Smith 1464b9ad928SBarry Smith #undef __FUNCT__ 1474b9ad928SBarry Smith #define __FUNCT__ "PCDestroy_Redundant" 148*6849ba73SBarry Smith static PetscErrorCode PCDestroy_Redundant(PC pc) 1494b9ad928SBarry Smith { 1504b9ad928SBarry Smith PC_Redundant *red = (PC_Redundant*)pc->data; 151dfbe8321SBarry Smith PetscErrorCode ierr; 1524b9ad928SBarry Smith 1534b9ad928SBarry Smith PetscFunctionBegin; 1544b9ad928SBarry Smith if (red->scatterin) {ierr = VecScatterDestroy(red->scatterin);CHKERRQ(ierr);} 1554b9ad928SBarry Smith if (red->scatterout) {ierr = VecScatterDestroy(red->scatterout);CHKERRQ(ierr);} 1564b9ad928SBarry Smith if (red->x) {ierr = VecDestroy(red->x);CHKERRQ(ierr);} 1574b9ad928SBarry Smith if (red->b) {ierr = VecDestroy(red->b);CHKERRQ(ierr);} 1584b9ad928SBarry Smith if (red->pmats) { 1594b9ad928SBarry Smith ierr = MatDestroyMatrices(1,&red->pmats);CHKERRQ(ierr); 1604b9ad928SBarry Smith } 1614b9ad928SBarry Smith ierr = PCDestroy(red->pc);CHKERRQ(ierr); 1624b9ad928SBarry Smith ierr = PetscFree(red);CHKERRQ(ierr); 1634b9ad928SBarry Smith PetscFunctionReturn(0); 1644b9ad928SBarry Smith } 1654b9ad928SBarry Smith 1664b9ad928SBarry Smith #undef __FUNCT__ 1674b9ad928SBarry Smith #define __FUNCT__ "PCSetFromOptions_Redundant" 168*6849ba73SBarry Smith static PetscErrorCode PCSetFromOptions_Redundant(PC pc) 1694b9ad928SBarry Smith { 1704b9ad928SBarry Smith PetscFunctionBegin; 1714b9ad928SBarry Smith PetscFunctionReturn(0); 1724b9ad928SBarry Smith } 1734b9ad928SBarry Smith 1744b9ad928SBarry Smith EXTERN_C_BEGIN 1754b9ad928SBarry Smith #undef __FUNCT__ 1764b9ad928SBarry Smith #define __FUNCT__ "PCRedundantSetScatter_Redundant" 177dfbe8321SBarry Smith PetscErrorCode PCRedundantSetScatter_Redundant(PC pc,VecScatter in,VecScatter out) 1784b9ad928SBarry Smith { 1794b9ad928SBarry Smith PC_Redundant *red = (PC_Redundant*)pc->data; 180dfbe8321SBarry Smith PetscErrorCode ierr; 1814b9ad928SBarry Smith 1824b9ad928SBarry Smith PetscFunctionBegin; 1834b9ad928SBarry Smith red->scatterin = in; 1844b9ad928SBarry Smith red->scatterout = out; 1854b9ad928SBarry Smith ierr = PetscObjectReference((PetscObject)in);CHKERRQ(ierr); 1864b9ad928SBarry Smith ierr = PetscObjectReference((PetscObject)out);CHKERRQ(ierr); 1874b9ad928SBarry Smith PetscFunctionReturn(0); 1884b9ad928SBarry Smith } 1894b9ad928SBarry Smith EXTERN_C_END 1904b9ad928SBarry Smith 1914b9ad928SBarry Smith #undef __FUNCT__ 1924b9ad928SBarry Smith #define __FUNCT__ "PCRedundantSetScatter" 1934b9ad928SBarry Smith /*@ 1944b9ad928SBarry Smith PCRedundantSetScatter - Sets the scatter used to copy values into the 1954b9ad928SBarry Smith redundant local solve and the scatter to move them back into the global 1964b9ad928SBarry Smith vector. 1974b9ad928SBarry Smith 1984b9ad928SBarry Smith Collective on PC 1994b9ad928SBarry Smith 2004b9ad928SBarry Smith Input Parameters: 2014b9ad928SBarry Smith + pc - the preconditioner context 2024b9ad928SBarry Smith . in - the scatter to move the values in 2034b9ad928SBarry Smith - out - the scatter to move them out 2044b9ad928SBarry Smith 2054b9ad928SBarry Smith Level: advanced 2064b9ad928SBarry Smith 2074b9ad928SBarry Smith .keywords: PC, redundant solve 2084b9ad928SBarry Smith @*/ 209dfbe8321SBarry Smith PetscErrorCode PCRedundantSetScatter(PC pc,VecScatter in,VecScatter out) 2104b9ad928SBarry Smith { 211dfbe8321SBarry Smith PetscErrorCode ierr,(*f)(PC,VecScatter,VecScatter); 2124b9ad928SBarry Smith 2134b9ad928SBarry Smith PetscFunctionBegin; 2144482741eSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 2154482741eSBarry Smith PetscValidHeaderSpecific(in,VEC_SCATTER_COOKIE,2); 2164482741eSBarry Smith PetscValidHeaderSpecific(out,VEC_SCATTER_COOKIE,3); 2174b9ad928SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCRedundantSetScatter_C",(void (**)(void))&f);CHKERRQ(ierr); 2184b9ad928SBarry Smith if (f) { 2194b9ad928SBarry Smith ierr = (*f)(pc,in,out);CHKERRQ(ierr); 2204b9ad928SBarry Smith } 2214b9ad928SBarry Smith PetscFunctionReturn(0); 2224b9ad928SBarry Smith } 2234b9ad928SBarry Smith 2244b9ad928SBarry Smith EXTERN_C_BEGIN 2254b9ad928SBarry Smith #undef __FUNCT__ 2264b9ad928SBarry Smith #define __FUNCT__ "PCRedundantGetPC_Redundant" 227dfbe8321SBarry Smith PetscErrorCode PCRedundantGetPC_Redundant(PC pc,PC *innerpc) 2284b9ad928SBarry Smith { 2294b9ad928SBarry Smith PC_Redundant *red = (PC_Redundant*)pc->data; 2304b9ad928SBarry Smith 2314b9ad928SBarry Smith PetscFunctionBegin; 2324b9ad928SBarry Smith *innerpc = red->pc; 2334b9ad928SBarry Smith PetscFunctionReturn(0); 2344b9ad928SBarry Smith } 2354b9ad928SBarry Smith EXTERN_C_END 2364b9ad928SBarry Smith 2374b9ad928SBarry Smith #undef __FUNCT__ 2384b9ad928SBarry Smith #define __FUNCT__ "PCRedundantGetPC" 2394b9ad928SBarry Smith /*@ 2404b9ad928SBarry Smith PCRedundantGetPC - Gets the sequential PC created by the redundant PC. 2414b9ad928SBarry Smith 2424b9ad928SBarry Smith Not Collective 2434b9ad928SBarry Smith 2444b9ad928SBarry Smith Input Parameter: 2454b9ad928SBarry Smith . pc - the preconditioner context 2464b9ad928SBarry Smith 2474b9ad928SBarry Smith Output Parameter: 2484b9ad928SBarry Smith . innerpc - the sequential PC 2494b9ad928SBarry Smith 2504b9ad928SBarry Smith Level: advanced 2514b9ad928SBarry Smith 2524b9ad928SBarry Smith .keywords: PC, redundant solve 2534b9ad928SBarry Smith @*/ 254dfbe8321SBarry Smith PetscErrorCode PCRedundantGetPC(PC pc,PC *innerpc) 2554b9ad928SBarry Smith { 256dfbe8321SBarry Smith PetscErrorCode ierr,(*f)(PC,PC*); 2574b9ad928SBarry Smith 2584b9ad928SBarry Smith PetscFunctionBegin; 2594482741eSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 2604482741eSBarry Smith PetscValidPointer(innerpc,2); 2614b9ad928SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCRedundantGetPC_C",(void (**)(void))&f);CHKERRQ(ierr); 2624b9ad928SBarry Smith if (f) { 2634b9ad928SBarry Smith ierr = (*f)(pc,innerpc);CHKERRQ(ierr); 2644b9ad928SBarry Smith } 2654b9ad928SBarry Smith PetscFunctionReturn(0); 2664b9ad928SBarry Smith } 2674b9ad928SBarry Smith 2684b9ad928SBarry Smith EXTERN_C_BEGIN 2694b9ad928SBarry Smith #undef __FUNCT__ 2704b9ad928SBarry Smith #define __FUNCT__ "PCRedundantGetOperators_Redundant" 271dfbe8321SBarry Smith PetscErrorCode PCRedundantGetOperators_Redundant(PC pc,Mat *mat,Mat *pmat) 2724b9ad928SBarry Smith { 2734b9ad928SBarry Smith PC_Redundant *red = (PC_Redundant*)pc->data; 2744b9ad928SBarry Smith 2754b9ad928SBarry Smith PetscFunctionBegin; 2764b9ad928SBarry Smith if (mat) *mat = red->pmats[0]; 2774b9ad928SBarry Smith if (pmat) *pmat = red->pmats[0]; 2784b9ad928SBarry Smith PetscFunctionReturn(0); 2794b9ad928SBarry Smith } 2804b9ad928SBarry Smith EXTERN_C_END 2814b9ad928SBarry Smith 2824b9ad928SBarry Smith #undef __FUNCT__ 2834b9ad928SBarry Smith #define __FUNCT__ "PCRedundantGetOperators" 2844b9ad928SBarry Smith /*@ 2854b9ad928SBarry Smith PCRedundantGetOperators - gets the sequential matrix and preconditioner matrix 2864b9ad928SBarry Smith 2874b9ad928SBarry Smith Not Collective 2884b9ad928SBarry Smith 2894b9ad928SBarry Smith Input Parameter: 2904b9ad928SBarry Smith . pc - the preconditioner context 2914b9ad928SBarry Smith 2924b9ad928SBarry Smith Output Parameters: 2934b9ad928SBarry Smith + mat - the matrix 2944b9ad928SBarry Smith - pmat - the (possibly different) preconditioner matrix 2954b9ad928SBarry Smith 2964b9ad928SBarry Smith Level: advanced 2974b9ad928SBarry Smith 2984b9ad928SBarry Smith .keywords: PC, redundant solve 2994b9ad928SBarry Smith @*/ 300dfbe8321SBarry Smith PetscErrorCode PCRedundantGetOperators(PC pc,Mat *mat,Mat *pmat) 3014b9ad928SBarry Smith { 302dfbe8321SBarry Smith PetscErrorCode ierr,(*f)(PC,Mat*,Mat*); 3034b9ad928SBarry Smith 3044b9ad928SBarry Smith PetscFunctionBegin; 3054482741eSBarry Smith PetscValidHeaderSpecific(pc,PC_COOKIE,1); 3064482741eSBarry Smith if (mat) PetscValidPointer(mat,2); 3074482741eSBarry Smith if (pmat) PetscValidPointer(pmat,3); 3084b9ad928SBarry Smith ierr = PetscObjectQueryFunction((PetscObject)pc,"PCRedundantGetOperators_C",(void (**)(void))&f);CHKERRQ(ierr); 3094b9ad928SBarry Smith if (f) { 3104b9ad928SBarry Smith ierr = (*f)(pc,mat,pmat);CHKERRQ(ierr); 3114b9ad928SBarry Smith } 3124b9ad928SBarry Smith PetscFunctionReturn(0); 3134b9ad928SBarry Smith } 3144b9ad928SBarry Smith 3154b9ad928SBarry Smith /* -------------------------------------------------------------------------------------*/ 31637a17b4dSBarry Smith /*MC 31737a17b4dSBarry Smith PCREDUNDANT - Runs a preconditioner for the entire problem on each processor 31837a17b4dSBarry Smith 31937a17b4dSBarry Smith 32037a17b4dSBarry Smith Options for the redundant preconditioners can be set with -redundant_pc_xxx 32137a17b4dSBarry Smith 32237a17b4dSBarry Smith Level: intermediate 32337a17b4dSBarry Smith 32437a17b4dSBarry Smith 32537a17b4dSBarry Smith .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PCRedundantSetScatter(), 32637a17b4dSBarry Smith PCRedundantGetPC(), PCRedundantGetOperators() 32737a17b4dSBarry Smith 32837a17b4dSBarry Smith M*/ 32937a17b4dSBarry Smith 3304b9ad928SBarry Smith EXTERN_C_BEGIN 3314b9ad928SBarry Smith #undef __FUNCT__ 3324b9ad928SBarry Smith #define __FUNCT__ "PCCreate_Redundant" 333dfbe8321SBarry Smith PetscErrorCode PCCreate_Redundant(PC pc) 3344b9ad928SBarry Smith { 335dfbe8321SBarry Smith PetscErrorCode ierr; 3364b9ad928SBarry Smith PC_Redundant *red; 3374b9ad928SBarry Smith char *prefix; 3384b9ad928SBarry Smith 3394b9ad928SBarry Smith PetscFunctionBegin; 3404b9ad928SBarry Smith ierr = PetscNew(PC_Redundant,&red);CHKERRQ(ierr); 3414b9ad928SBarry Smith PetscLogObjectMemory(pc,sizeof(PC_Redundant)); 3424b9ad928SBarry Smith ierr = PetscMemzero(red,sizeof(PC_Redundant));CHKERRQ(ierr); 3434b9ad928SBarry Smith red->useparallelmat = PETSC_TRUE; 3444b9ad928SBarry Smith 3454b9ad928SBarry Smith /* create the sequential PC that each processor has copy of */ 3464b9ad928SBarry Smith ierr = PCCreate(PETSC_COMM_SELF,&red->pc);CHKERRQ(ierr); 3474b9ad928SBarry Smith ierr = PCSetType(red->pc,PCLU);CHKERRQ(ierr); 3484b9ad928SBarry Smith ierr = PCGetOptionsPrefix(pc,&prefix);CHKERRQ(ierr); 3494b9ad928SBarry Smith ierr = PCSetOptionsPrefix(red->pc,prefix);CHKERRQ(ierr); 3504b9ad928SBarry Smith ierr = PCAppendOptionsPrefix(red->pc,"redundant_");CHKERRQ(ierr); 3514b9ad928SBarry Smith 3524b9ad928SBarry Smith pc->ops->apply = PCApply_Redundant; 3534b9ad928SBarry Smith pc->ops->applytranspose = 0; 3544b9ad928SBarry Smith pc->ops->setup = PCSetUp_Redundant; 3554b9ad928SBarry Smith pc->ops->destroy = PCDestroy_Redundant; 3564b9ad928SBarry Smith pc->ops->setfromoptions = PCSetFromOptions_Redundant; 3574b9ad928SBarry Smith pc->ops->setuponblocks = 0; 3584b9ad928SBarry Smith pc->ops->view = PCView_Redundant; 3594b9ad928SBarry Smith pc->ops->applyrichardson = 0; 3604b9ad928SBarry Smith 3614b9ad928SBarry Smith pc->data = (void*)red; 3624b9ad928SBarry Smith 3634b9ad928SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCRedundantSetScatter_C","PCRedundantSetScatter_Redundant", 3644b9ad928SBarry Smith PCRedundantSetScatter_Redundant);CHKERRQ(ierr); 3654b9ad928SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCRedundantGetPC_C","PCRedundantGetPC_Redundant", 3664b9ad928SBarry Smith PCRedundantGetPC_Redundant);CHKERRQ(ierr); 3674b9ad928SBarry Smith ierr = PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCRedundantGetOperators_C","PCRedundantGetOperators_Redundant", 3684b9ad928SBarry Smith PCRedundantGetOperators_Redundant);CHKERRQ(ierr); 3694b9ad928SBarry Smith 3704b9ad928SBarry Smith PetscFunctionReturn(0); 3714b9ad928SBarry Smith } 3724b9ad928SBarry Smith EXTERN_C_END 373