15b89ad90SMark F. Adams /* 25b89ad90SMark F. Adams GAMG geometric-algebric multiogrid PC - Mark Adams 2011 35b89ad90SMark F. Adams */ 45b89ad90SMark F. Adams #include <private/pcimpl.h> /*I "petscpc.h" I*/ 55b89ad90SMark F. Adams #include <../src/ksp/pc/impls/mg/mgimpl.h> /*I "petscpcmg.h" I*/ 65b89ad90SMark F. Adams #include <../src/mat/impls/aij/seq/aij.h> 75b89ad90SMark F. Adams #include <../src/mat/impls/aij/mpi/mpiaij.h> 85b89ad90SMark F. Adams 95b89ad90SMark F. Adams /* Private context for the GAMG preconditioner */ 105b89ad90SMark F. Adams typedef struct gamg_TAG { 115b89ad90SMark F. Adams PetscInt m_dim; 125b89ad90SMark F. Adams PetscInt m_Nlevels; 135b89ad90SMark F. Adams PetscInt m_data_sz; 145b89ad90SMark F. Adams PetscReal *m_data; /* blocked vector of vertex data on fine grid (coordinates) */ 155b89ad90SMark F. Adams } PC_GAMG; 165b89ad90SMark F. Adams 175b89ad90SMark F. Adams /* -----------------------------------------------------------------------------*/ 185b89ad90SMark F. Adams #undef __FUNCT__ 195b89ad90SMark F. Adams #define __FUNCT__ "PCReset_GAMG" 205b89ad90SMark F. Adams PetscErrorCode PCReset_GAMG(PC pc) 215b89ad90SMark F. Adams { 225b89ad90SMark F. Adams PetscErrorCode ierr; 235b89ad90SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 245b89ad90SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 255b89ad90SMark F. Adams 265b89ad90SMark F. Adams PetscFunctionBegin; 275b89ad90SMark F. Adams ierr = PetscFree(pc_gamg->m_data);CHKERRQ(ierr); 285b89ad90SMark F. Adams PetscFunctionReturn(0); 295b89ad90SMark F. Adams } 305b89ad90SMark F. Adams 315b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */ 325b89ad90SMark F. Adams /* 335b89ad90SMark F. Adams PCSetCoordinates_GAMG 345b89ad90SMark F. Adams 355b89ad90SMark F. Adams Input Parameter: 365b89ad90SMark F. Adams . pc - the preconditioner context 375b89ad90SMark F. Adams */ 385b89ad90SMark F. Adams #undef __FUNCT__ 395b89ad90SMark F. Adams #define __FUNCT__ "PCSetCoordinates_GAMG" 40*6c237d78SBarry Smith PetscErrorCode PCSetCoordinates_GAMG(PC pc, const int ndm,PetscReal *coords ) 415b89ad90SMark F. Adams { 425b89ad90SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 435b89ad90SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 44*6c237d78SBarry Smith PetscErrorCode ierr; 45*6c237d78SBarry Smith PetscInt bs, my0, tt; 46*6c237d78SBarry Smith Mat mat = pc->pmat; 47*6c237d78SBarry Smith PetscInt arrsz; 485b89ad90SMark F. Adams 495b89ad90SMark F. Adams PetscFunctionBegin; 505b89ad90SMark F. Adams ierr = MatGetBlockSize( mat, &bs ); CHKERRQ( ierr ); 515b89ad90SMark F. Adams ierr = MatGetOwnershipRange( mat, &my0, &tt ); CHKERRQ(ierr); 52*6c237d78SBarry Smith arrsz = (tt-my0)/bs*ndm; 535b89ad90SMark F. Adams 545b89ad90SMark F. Adams // put coordinates 55*6c237d78SBarry Smith if (!pc_gamg->m_data || (pc_gamg->m_data_sz != arrsz)) { 565b89ad90SMark F. Adams ierr = PetscFree( pc_gamg->m_data ); CHKERRQ(ierr); 57*6c237d78SBarry Smith ierr = PetscMalloc(arrsz*sizeof(double), &pc_gamg->m_data ); CHKERRQ(ierr); 585b89ad90SMark F. Adams } 595b89ad90SMark F. Adams 605b89ad90SMark F. Adams /* copy data in */ 615b89ad90SMark F. Adams for(tt=0;tt<arrsz;tt++){ 625b89ad90SMark F. Adams pc_gamg->m_data[tt] = coords[tt]; 635b89ad90SMark F. Adams } 645b89ad90SMark F. Adams pc_gamg->m_data_sz = arrsz; 655b89ad90SMark F. Adams pc_gamg->m_dim = ndm; 665b89ad90SMark F. Adams 675b89ad90SMark F. Adams PetscFunctionReturn(0); 685b89ad90SMark F. Adams } 695b89ad90SMark F. Adams 705b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */ 715b89ad90SMark F. Adams /* 725b89ad90SMark F. Adams createCrsOp 735b89ad90SMark F. Adams 745b89ad90SMark F. Adams Input Parameter: 755b89ad90SMark F. Adams . Amat - matrix on this fine level 765b89ad90SMark F. Adams . P_out - prolongation operator to the next level 775b89ad90SMark F. Adams . Acrs - coarse matrix that is created 785b89ad90SMark F. Adams */ 795b89ad90SMark F. Adams #undef __FUNCT__ 805b89ad90SMark F. Adams #define __FUNCT__ "createCrsOp" 815b89ad90SMark F. Adams PetscErrorCode createCrsOp( Mat Amat, Mat P_inout, Mat *Acrs ) 825b89ad90SMark F. Adams { 835b89ad90SMark F. Adams PetscErrorCode ierr; 845b89ad90SMark F. Adams 855b89ad90SMark F. Adams PetscFunctionBegin; 865b89ad90SMark F. Adams Mat H; 875b89ad90SMark F. Adams ierr = MatPtAP( Amat, P_inout, MAT_INITIAL_MATRIX, 2.0, &H); CHKERRQ(ierr); 88acadaa72SMark F. Adams 89acadaa72SMark F. Adams /* need to repartition H and move colums of P accordingly */ 90acadaa72SMark F. Adams 91acadaa72SMark F. Adams 925b89ad90SMark F. Adams *Acrs = H; 935b89ad90SMark F. Adams 945b89ad90SMark F. Adams PetscFunctionReturn(0); 955b89ad90SMark F. Adams } 965b89ad90SMark F. Adams 975b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */ 985b89ad90SMark F. Adams /* 995b89ad90SMark F. Adams PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner 1005b89ad90SMark F. Adams by setting data structures and options. 1015b89ad90SMark F. Adams 1025b89ad90SMark F. Adams Input Parameter: 1035b89ad90SMark F. Adams . pc - the preconditioner context 1045b89ad90SMark F. Adams 1055b89ad90SMark F. Adams Application Interface Routine: PCSetUp() 1065b89ad90SMark F. Adams 1075b89ad90SMark F. Adams Notes: 1085b89ad90SMark F. Adams The interface routine PCSetUp() is not usually called directly by 1095b89ad90SMark F. Adams the user, but instead is called by PCApply() if necessary. 1105b89ad90SMark F. Adams */ 1115b89ad90SMark F. Adams extern PetscErrorCode PCSetFromOptions_MG(PC); 1125b89ad90SMark F. Adams extern PetscErrorCode PCReset_MG(PC); 1135b89ad90SMark F. Adams extern PetscErrorCode createProlongation( Mat, Mat *, PetscReal [], PetscReal **, const PetscInt); 1145b89ad90SMark F. Adams #undef __FUNCT__ 1155b89ad90SMark F. Adams #define __FUNCT__ "PCSetUp_GAMG" 1165b89ad90SMark F. Adams PetscErrorCode PCSetUp_GAMG( PC pc ) 1175b89ad90SMark F. Adams { 1185b89ad90SMark F. Adams PetscErrorCode ierr; 1195b89ad90SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 1205b89ad90SMark F. Adams PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 1215b89ad90SMark F. Adams Mat Amat = pc->mat, Pmat = pc->pmat; 1225b89ad90SMark F. Adams PetscBool isSeq, isMPI; 1235b89ad90SMark F. Adams PetscInt fine_level, level, level1, M, N, bs, lidx; 1245b89ad90SMark F. Adams MPI_Comm wcomm = ((PetscObject)pc)->comm; 1255b89ad90SMark F. Adams 1265b89ad90SMark F. Adams PetscFunctionBegin; 1275b89ad90SMark F. Adams if (pc->setupcalled){ 1285b89ad90SMark F. Adams /* no state data in GAMG to destroy (now) */ 1295b89ad90SMark F. Adams ierr = PCReset_MG(pc);CHKERRQ(ierr); 1305b89ad90SMark F. Adams } 131*6c237d78SBarry Smith if (!pc_gamg->m_data) SETERRQ(wcomm,PETSC_ERR_SUP,"PCSetUp_GAMG called before PCSetCoordinates"); 1325b89ad90SMark F. Adams /* setup special features of PCGAMG */ 1335b89ad90SMark F. Adams ierr = PetscTypeCompare((PetscObject) Amat, MATSEQAIJ, &isSeq);CHKERRQ(ierr); 1345b89ad90SMark F. Adams ierr = PetscTypeCompare((PetscObject) Amat, MATMPIAIJ, &isMPI);CHKERRQ(ierr); 1355b89ad90SMark F. Adams if (isMPI) { 1365b89ad90SMark F. Adams } else if (isSeq) { 1375b89ad90SMark F. Adams } else SETERRQ1(wcomm,PETSC_ERR_ARG_WRONG, "Matrix type '%s' cannot be used with GAMG. GAMG can only handle AIJ matrices.",((PetscObject)Amat)->type_name); 1385b89ad90SMark F. Adams 1395b89ad90SMark F. Adams /* GAMG requires input of fine-grid matrix. It determines nlevels. */ 1405b89ad90SMark F. Adams ierr = MatGetSize( Amat, &M, &N );CHKERRQ(ierr); 1415b89ad90SMark F. Adams ierr = MatGetBlockSize( Amat, &bs ); CHKERRQ(ierr); 1425b89ad90SMark F. Adams if(bs!=1) SETERRQ1(wcomm,PETSC_ERR_ARG_WRONG, "GAMG only supports scalar prblems bs = '%d'.",bs); 1435b89ad90SMark F. Adams 1445b89ad90SMark F. Adams /* Get A_i and R_i */ 1455b89ad90SMark F. Adams #define GAMG_MAXLEVELS 10 1465b89ad90SMark F. Adams Mat Aarr[GAMG_MAXLEVELS], Rarr[GAMG_MAXLEVELS]; PetscReal *coarse_crds = 0, *crds = pc_gamg->m_data; 1475b89ad90SMark F. Adams for (level=0, Aarr[0] = Pmat; level < GAMG_MAXLEVELS-1; level++ ){ 1485b89ad90SMark F. Adams ierr = MatGetSize( Aarr[level], &M, &N );CHKERRQ(ierr); 1495b89ad90SMark F. Adams if( M < 100 ) { /* hard wire this for now */ 1505b89ad90SMark F. Adams break; 1515b89ad90SMark F. Adams } 1525b89ad90SMark F. Adams level1 = level + 1; 1535b89ad90SMark F. Adams ierr = createProlongation( Aarr[level], &Rarr[level1], crds, &coarse_crds, pc_gamg->m_dim ); 1545b89ad90SMark F. Adams CHKERRQ(ierr); 1555b89ad90SMark F. Adams if(level==0) Aarr[0] = Amat; /* use Pmat for finest level setup, but use mat for solver */ 1565b89ad90SMark F. Adams ierr = createCrsOp( Aarr[level], Rarr[level1], &Aarr[level1] ); CHKERRQ(ierr); 1575b89ad90SMark F. Adams ierr = PetscFree( crds ); CHKERRQ( ierr ); 1585b89ad90SMark F. Adams crds = coarse_crds; 1595b89ad90SMark F. Adams } 1605b89ad90SMark F. Adams ierr = PetscFree( coarse_crds ); CHKERRQ( ierr ); 161*6c237d78SBarry Smith 1625b89ad90SMark F. Adams pc_gamg->m_data = 0; /* destroyed coordinate data */ 1635b89ad90SMark F. Adams pc_gamg->m_Nlevels = level + 1; 1645b89ad90SMark F. Adams fine_level = level; 1655b89ad90SMark F. Adams ierr = PCMGSetLevels(pc,pc_gamg->m_Nlevels,PETSC_NULL);CHKERRQ(ierr); 1665b89ad90SMark F. Adams 1675b89ad90SMark F. Adams /* set default smoothers */ 1685b89ad90SMark F. Adams PetscReal emax = 2.0, emin; 1695b89ad90SMark F. Adams for (level=1; level<=fine_level; level++){ 1705b89ad90SMark F. Adams KSP smoother; PC subpc; 1715b89ad90SMark F. Adams ierr = PCMGGetSmoother(pc,level,&smoother);CHKERRQ(ierr); 1725b89ad90SMark F. Adams ierr = KSPSetType(smoother,KSPCHEBYCHEV);CHKERRQ(ierr); 1735b89ad90SMark F. Adams emin = emax/5.0; 1745b89ad90SMark F. Adams ierr = KSPGetPC(smoother,&subpc);CHKERRQ(ierr); 1755b89ad90SMark F. Adams ierr = PCSetType(subpc,PCJACOBI);CHKERRQ(ierr); 1765b89ad90SMark F. Adams ierr = KSPChebychevSetEigenvalues(smoother, emax, emin);CHKERRQ(ierr); /* need auto !!!!*/ 1775b89ad90SMark F. Adams } 1785b89ad90SMark F. Adams ierr = PCSetFromOptions_MG(pc); CHKERRQ(ierr); /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */ 1795b89ad90SMark F. Adams { 1805b89ad90SMark F. Adams PetscBool galerkin; 1815b89ad90SMark F. Adams ierr = PCMGGetGalerkin( pc, &galerkin); CHKERRQ(ierr); 1825b89ad90SMark F. Adams if(galerkin){ 1835b89ad90SMark F. Adams SETERRQ(wcomm,PETSC_ERR_ARG_WRONG, "GAMG does galerkin manually so it must not be used in PC_MG."); 1845b89ad90SMark F. Adams } 1855b89ad90SMark F. Adams } 1865b89ad90SMark F. Adams /* create coarse level and the interpolation between the levels */ 1875b89ad90SMark F. Adams for (level=0,lidx=pc_gamg->m_Nlevels-1; level<fine_level; level++,lidx--){ 1885b89ad90SMark F. Adams level1 = level + 1; 1895b89ad90SMark F. Adams /* PetscInt MM,NN; */ 1905b89ad90SMark F. Adams /* ierr = MatGetSize( Rarr[lidx], &MM, &NN );CHKERRQ(ierr); */ 1915b89ad90SMark F. Adams /* PetscPrintf(PETSC_COMM_WORLD,"%s Set P(%d,%d) on level %d (%d)\n",__FUNCT__,MM,NN,level1,lidx); */ 1925b89ad90SMark F. Adams ierr = PCMGSetInterpolation(pc,level1,Rarr[lidx]);CHKERRQ(ierr); 193*6c237d78SBarry Smith if(!PETSC_TRUE) { 1945b89ad90SMark F. Adams PetscViewer viewer; 1955b89ad90SMark F. Adams ierr = PetscViewerASCIIOpen(PETSC_COMM_SELF, "Rmat.m", &viewer); CHKERRQ(ierr); 1965b89ad90SMark F. Adams ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB); CHKERRQ(ierr); 1975b89ad90SMark F. Adams ierr = MatView(Rarr[level1],viewer);CHKERRQ(ierr); 1985b89ad90SMark F. Adams ierr = PetscViewerDestroy( &viewer ); 1995b89ad90SMark F. Adams } 2005b89ad90SMark F. Adams KSP smoother; 2015b89ad90SMark F. Adams ierr = PCMGGetSmoother(pc,level,&smoother); CHKERRQ(ierr); 2025b89ad90SMark F. Adams ierr = KSPSetOperators( smoother, Aarr[lidx], Aarr[lidx], DIFFERENT_NONZERO_PATTERN ); 2035b89ad90SMark F. Adams CHKERRQ(ierr); 2045b89ad90SMark F. Adams /* ierr = MatGetSize( Aarr[lidx], &MM, &NN ); CHKERRQ(ierr); */ 2055b89ad90SMark F. Adams /* PetscPrintf(PETSC_COMM_WORLD,"%s Set A(%d,%d) on level %d (%d)\n",__FUNCT__,MM,NN,level,lidx); */ 2065b89ad90SMark F. Adams } 2075b89ad90SMark F. Adams { /* fine level (no P) */ 2085b89ad90SMark F. Adams KSP smoother; 2095b89ad90SMark F. Adams ierr = PCMGGetSmoother(pc,fine_level,&smoother); CHKERRQ(ierr); 2105b89ad90SMark F. Adams ierr = KSPSetOperators( smoother, Aarr[0], Aarr[0], DIFFERENT_NONZERO_PATTERN ); 2115b89ad90SMark F. Adams CHKERRQ(ierr); 2125b89ad90SMark F. Adams /* PetscInt MM,NN; */ 2135b89ad90SMark F. Adams /* ierr = MatGetSize( Aarr[0], &MM, &NN );CHKERRQ(ierr); */ 2145b89ad90SMark F. Adams /* PetscPrintf(PETSC_COMM_WORLD,"%s Set A(%d,%d) on level %d (%d)\n",__FUNCT__,MM,NN,fine_level,0); */ 2155b89ad90SMark F. Adams } 2165b89ad90SMark F. Adams 2175b89ad90SMark F. Adams /* setupcalled is set to 0 so that MG is setup from scratch */ 2185b89ad90SMark F. Adams pc->setupcalled = 0; 2195b89ad90SMark F. Adams ierr = PCSetUp_MG(pc);CHKERRQ(ierr); 2205b89ad90SMark F. Adams PetscFunctionReturn(0); 2215b89ad90SMark F. Adams } 2225b89ad90SMark F. Adams 2235b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */ 2245b89ad90SMark F. Adams /* 2255b89ad90SMark F. Adams PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner 2265b89ad90SMark F. Adams that was created with PCCreate_GAMG(). 2275b89ad90SMark F. Adams 2285b89ad90SMark F. Adams Input Parameter: 2295b89ad90SMark F. Adams . pc - the preconditioner context 2305b89ad90SMark F. Adams 2315b89ad90SMark F. Adams Application Interface Routine: PCDestroy() 2325b89ad90SMark F. Adams */ 2335b89ad90SMark F. Adams #undef __FUNCT__ 2345b89ad90SMark F. Adams #define __FUNCT__ "PCDestroy_GAMG" 2355b89ad90SMark F. Adams PetscErrorCode PCDestroy_GAMG(PC pc) 2365b89ad90SMark F. Adams { 2375b89ad90SMark F. Adams PetscErrorCode ierr; 2385b89ad90SMark F. Adams PC_MG *mg = (PC_MG*)pc->data; 2395b89ad90SMark F. Adams PC_GAMG *pc_gamg= (PC_GAMG*)mg->innerctx; 2405b89ad90SMark F. Adams 2415b89ad90SMark F. Adams PetscFunctionBegin; 2425b89ad90SMark F. Adams ierr = PCReset_GAMG(pc);CHKERRQ(ierr); 2435b89ad90SMark F. Adams ierr = PetscFree(pc_gamg);CHKERRQ(ierr); 2445b89ad90SMark F. Adams ierr = PCDestroy_MG(pc);CHKERRQ(ierr); 2455b89ad90SMark F. Adams PetscFunctionReturn(0); 2465b89ad90SMark F. Adams } 2475b89ad90SMark F. Adams 2485b89ad90SMark F. Adams #undef __FUNCT__ 2495b89ad90SMark F. Adams #define __FUNCT__ "PCSetFromOptions_GAMG" 2505b89ad90SMark F. Adams PetscErrorCode PCSetFromOptions_GAMG(PC pc) 2515b89ad90SMark F. Adams { 2525b89ad90SMark F. Adams /* PetscErrorCode ierr; */ 2535b89ad90SMark F. Adams /* PC_MG *mg = (PC_MG*)pc->data; */ 2545b89ad90SMark F. Adams /* PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; */ 2555b89ad90SMark F. Adams /* MPI_Comm comm = ((PetscObject)pc)->comm; */ 2565b89ad90SMark F. Adams 2575b89ad90SMark F. Adams PetscFunctionBegin; 2585b89ad90SMark F. Adams PetscFunctionReturn(0); 2595b89ad90SMark F. Adams } 2605b89ad90SMark F. Adams 2615b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */ 2625b89ad90SMark F. Adams /* 2635b89ad90SMark F. Adams PCCreate_GAMG - Creates a GAMG preconditioner context, PC_GAMG 2645b89ad90SMark F. Adams 2655b89ad90SMark F. Adams Input Parameter: 2665b89ad90SMark F. Adams . pc - the preconditioner context 2675b89ad90SMark F. Adams 2685b89ad90SMark F. Adams Application Interface Routine: PCCreate() 2695b89ad90SMark F. Adams 2705b89ad90SMark F. Adams */ 2715b89ad90SMark F. Adams /* MC 2725b89ad90SMark F. Adams PCGAMG - Use algebraic multigrid preconditioning. This preconditioner requires you provide 2735b89ad90SMark F. Adams fine grid discretization matrix and coordinates on the fine grid. 2745b89ad90SMark F. Adams 2755b89ad90SMark F. Adams Options Database Key: 2765b89ad90SMark F. Adams Multigrid options(inherited) 2775b89ad90SMark F. Adams + -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (MGSetCycles) 2785b89ad90SMark F. Adams . -pc_mg_smoothup <1>: Number of post-smoothing steps (MGSetNumberSmoothUp) 2795b89ad90SMark F. Adams . -pc_mg_smoothdown <1>: Number of pre-smoothing steps (MGSetNumberSmoothDown) 2805b89ad90SMark F. Adams -pc_mg_type <multiplicative>: (one of) additive multiplicative full cascade kascade 2815b89ad90SMark F. Adams GAMG options: 2825b89ad90SMark F. Adams 2835b89ad90SMark F. Adams Level: intermediate 2845b89ad90SMark F. Adams Concepts: multigrid 2855b89ad90SMark F. Adams 2865b89ad90SMark F. Adams .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, 2875b89ad90SMark F. Adams PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), MPSetCycles(), PCMGSetNumberSmoothDown(), 2885b89ad90SMark F. Adams PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(), 2895b89ad90SMark F. Adams PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(), 2905b89ad90SMark F. Adams PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR() 2915b89ad90SMark F. Adams M */ 2925b89ad90SMark F. Adams 2935b89ad90SMark F. Adams EXTERN_C_BEGIN 2945b89ad90SMark F. Adams #undef __FUNCT__ 2955b89ad90SMark F. Adams #define __FUNCT__ "PCCreate_GAMG" 2965b89ad90SMark F. Adams PetscErrorCode PCCreate_GAMG(PC pc) 2975b89ad90SMark F. Adams { 2985b89ad90SMark F. Adams PetscErrorCode ierr; 2995b89ad90SMark F. Adams PC_GAMG *pc_gamg; 3005b89ad90SMark F. Adams PC_MG *mg; 3015b89ad90SMark F. Adams 3025b89ad90SMark F. Adams PetscFunctionBegin; 3035b89ad90SMark F. Adams /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */ 3045b89ad90SMark F. Adams ierr = PCSetType(pc,PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */ 3055b89ad90SMark F. Adams ierr = PetscObjectChangeTypeName((PetscObject)pc,PCGAMG);CHKERRQ(ierr); 3065b89ad90SMark F. Adams 3075b89ad90SMark F. Adams /* create a supporting struct and attach it to pc */ 3085b89ad90SMark F. Adams ierr = PetscNewLog(pc,PC_GAMG,&pc_gamg);CHKERRQ(ierr); 3095b89ad90SMark F. Adams mg = (PC_MG*)pc->data; 3105b89ad90SMark F. Adams mg->innerctx = pc_gamg; 3115b89ad90SMark F. Adams 3125b89ad90SMark F. Adams pc_gamg->m_Nlevels = -1; 3135b89ad90SMark F. Adams 3145b89ad90SMark F. Adams /* overwrite the pointers of PCMG by the functions of PCGAMG */ 3155b89ad90SMark F. Adams pc->ops->setfromoptions = PCSetFromOptions_GAMG; 3165b89ad90SMark F. Adams pc->ops->setup = PCSetUp_GAMG; 3175b89ad90SMark F. Adams pc->ops->reset = PCReset_GAMG; 3185b89ad90SMark F. Adams pc->ops->destroy = PCDestroy_GAMG; 3195b89ad90SMark F. Adams 3205b89ad90SMark F. Adams ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc, 3215b89ad90SMark F. Adams "PCSetCoordinates_C", 3225b89ad90SMark F. Adams "PCSetCoordinates_GAMG", 3235b89ad90SMark F. Adams PCSetCoordinates_GAMG);CHKERRQ(ierr); 3245b89ad90SMark F. Adams PetscFunctionReturn(0); 3255b89ad90SMark F. Adams } 3265b89ad90SMark F. Adams EXTERN_C_END 327