xref: /petsc/src/ksp/pc/impls/gamg/gamg.c (revision baa4e9fa94decc45528f349d244fecd9d9d631ec)
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 */
38a92563c5SMark F. Adams EXTERN_C_BEGIN
395b89ad90SMark F. Adams #undef __FUNCT__
405b89ad90SMark F. Adams #define __FUNCT__ "PCSetCoordinates_GAMG"
416c237d78SBarry Smith PetscErrorCode PCSetCoordinates_GAMG(PC pc, const int ndm,PetscReal *coords )
425b89ad90SMark F. Adams {
435b89ad90SMark F. Adams   PC_MG          *mg = (PC_MG*)pc->data;
445b89ad90SMark F. Adams   PC_GAMG        *pc_gamg = (PC_GAMG*)mg->innerctx;
456c237d78SBarry Smith   PetscErrorCode ierr;
466c237d78SBarry Smith   PetscInt       bs, my0, tt;
476c237d78SBarry Smith   Mat            mat = pc->pmat;
486c237d78SBarry Smith   PetscInt       arrsz;
495b89ad90SMark F. Adams 
505b89ad90SMark F. Adams   PetscFunctionBegin;
515b89ad90SMark F. Adams   ierr  = MatGetBlockSize( mat, &bs );               CHKERRQ( ierr );
525b89ad90SMark F. Adams   ierr  = MatGetOwnershipRange( mat, &my0, &tt ); CHKERRQ(ierr);
536c237d78SBarry Smith   arrsz = (tt-my0)/bs*ndm;
545b89ad90SMark F. Adams 
555b89ad90SMark F. Adams   // put coordinates
566c237d78SBarry Smith   if (!pc_gamg->m_data || (pc_gamg->m_data_sz != arrsz)) {
575b89ad90SMark F. Adams     ierr = PetscFree( pc_gamg->m_data );  CHKERRQ(ierr);
586c237d78SBarry Smith     ierr = PetscMalloc(arrsz*sizeof(double), &pc_gamg->m_data ); CHKERRQ(ierr);
595b89ad90SMark F. Adams   }
605b89ad90SMark F. Adams 
615b89ad90SMark F. Adams   /* copy data in */
625b89ad90SMark F. Adams   for(tt=0;tt<arrsz;tt++){
635b89ad90SMark F. Adams     pc_gamg->m_data[tt] = coords[tt];
645b89ad90SMark F. Adams   }
655b89ad90SMark F. Adams   pc_gamg->m_data_sz = arrsz;
665b89ad90SMark F. Adams   pc_gamg->m_dim = ndm;
675b89ad90SMark F. Adams 
685b89ad90SMark F. Adams   PetscFunctionReturn(0);
695b89ad90SMark F. Adams }
70a92563c5SMark F. Adams EXTERN_C_END
715b89ad90SMark F. Adams 
725b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
735b89ad90SMark F. Adams /*
745b89ad90SMark F. Adams    createCrsOp
755b89ad90SMark F. Adams 
765b89ad90SMark F. Adams    Input Parameter:
775b89ad90SMark F. Adams    . Amat - matrix on this fine level
785b89ad90SMark F. Adams    . P_out - prolongation operator to the next level
795b89ad90SMark F. Adams    . Acrs - coarse matrix that is created
805b89ad90SMark F. Adams */
815b89ad90SMark F. Adams #undef __FUNCT__
825b89ad90SMark F. Adams #define __FUNCT__ "createCrsOp"
835b89ad90SMark F. Adams PetscErrorCode createCrsOp( Mat Amat, Mat P_inout, Mat *Acrs )
845b89ad90SMark F. Adams {
855b89ad90SMark F. Adams   PetscErrorCode ierr;
865b89ad90SMark F. Adams 
875b89ad90SMark F. Adams   PetscFunctionBegin;
885b89ad90SMark F. Adams   Mat H;
895b89ad90SMark F. Adams   ierr = MatPtAP( Amat, P_inout, MAT_INITIAL_MATRIX, 2.0, &H); CHKERRQ(ierr);
90acadaa72SMark F. Adams 
91acadaa72SMark F. Adams   /* need to repartition H and move colums of P accordingly */
92acadaa72SMark F. Adams 
93acadaa72SMark F. Adams 
945b89ad90SMark F. Adams   *Acrs = H;
955b89ad90SMark F. Adams 
965b89ad90SMark F. Adams   PetscFunctionReturn(0);
975b89ad90SMark F. Adams }
985b89ad90SMark F. Adams 
995b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
1005b89ad90SMark F. Adams /*
1015b89ad90SMark F. Adams    PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner
1025b89ad90SMark F. Adams                     by setting data structures and options.
1035b89ad90SMark F. Adams 
1045b89ad90SMark F. Adams    Input Parameter:
1055b89ad90SMark F. Adams .  pc - the preconditioner context
1065b89ad90SMark F. Adams 
1075b89ad90SMark F. Adams    Application Interface Routine: PCSetUp()
1085b89ad90SMark F. Adams 
1095b89ad90SMark F. Adams    Notes:
1105b89ad90SMark F. Adams    The interface routine PCSetUp() is not usually called directly by
1115b89ad90SMark F. Adams    the user, but instead is called by PCApply() if necessary.
1125b89ad90SMark F. Adams */
1135b89ad90SMark F. Adams extern PetscErrorCode PCSetFromOptions_MG(PC);
1145b89ad90SMark F. Adams extern PetscErrorCode PCReset_MG(PC);
1150f9369f8SMark F. Adams extern PetscErrorCode createProlongation( Mat, PetscReal [], const PetscInt,
116*baa4e9faSMark F. Adams                                           Mat *, PetscReal **, PetscBool *a_isOK );
1175b89ad90SMark F. Adams #undef __FUNCT__
1185b89ad90SMark F. Adams #define __FUNCT__ "PCSetUp_GAMG"
1195b89ad90SMark F. Adams PetscErrorCode PCSetUp_GAMG( PC pc )
1205b89ad90SMark F. Adams {
1215b89ad90SMark F. Adams   PetscErrorCode  ierr;
1225b89ad90SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
1235b89ad90SMark F. Adams   PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx;
1245b89ad90SMark F. Adams   Mat              Amat = pc->mat, Pmat = pc->pmat;
1255b89ad90SMark F. Adams   PetscBool        isSeq, isMPI;
1265b89ad90SMark F. Adams   PetscInt         fine_level, level, level1, M, N, bs, lidx;
1275b89ad90SMark F. Adams   MPI_Comm         wcomm = ((PetscObject)pc)->comm;
128*baa4e9faSMark F. Adams   PetscMPIInt      mype,npe;
1295b89ad90SMark F. Adams 
1305b89ad90SMark F. Adams   PetscFunctionBegin;
131*baa4e9faSMark F. Adams   ierr = MPI_Comm_rank(wcomm,&mype);CHKERRQ(ierr);
132*baa4e9faSMark F. Adams   ierr = MPI_Comm_size(wcomm,&npe);CHKERRQ(ierr);
1335b89ad90SMark F. Adams   if (pc->setupcalled){
1345b89ad90SMark F. Adams     /* no state data in GAMG to destroy (now) */
1355b89ad90SMark F. Adams     ierr = PCReset_MG(pc); CHKERRQ(ierr);
1365b89ad90SMark F. Adams   }
1376c237d78SBarry Smith   if (!pc_gamg->m_data) SETERRQ(wcomm,PETSC_ERR_SUP,"PCSetUp_GAMG called before PCSetCoordinates");
1385b89ad90SMark F. Adams   /* setup special features of PCGAMG */
1395b89ad90SMark F. Adams   ierr = PetscTypeCompare((PetscObject)Amat, MATSEQAIJ, &isSeq);CHKERRQ(ierr);
1405b89ad90SMark F. Adams   ierr = PetscTypeCompare((PetscObject)Amat, MATMPIAIJ, &isMPI);CHKERRQ(ierr);
1415b89ad90SMark F. Adams   if (isMPI) {
1425b89ad90SMark F. Adams   } else if (isSeq) {
1435b89ad90SMark 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);
1445b89ad90SMark F. Adams 
1455b89ad90SMark F. Adams   /* GAMG requires input of fine-grid matrix. It determines nlevels. */
1465b89ad90SMark F. Adams   ierr = MatGetSize( Amat, &M, &N );CHKERRQ(ierr);
1475b89ad90SMark F. Adams   ierr = MatGetBlockSize( Amat, &bs ); CHKERRQ(ierr);
1485b89ad90SMark F. Adams   if(bs!=1) SETERRQ1(wcomm,PETSC_ERR_ARG_WRONG, "GAMG only supports scalar prblems bs = '%d'.",bs);
1495b89ad90SMark F. Adams 
1505b89ad90SMark F. Adams   /* Get A_i and R_i */
1515b89ad90SMark F. Adams #define GAMG_MAXLEVELS 10
1525b89ad90SMark F. Adams   Mat Aarr[GAMG_MAXLEVELS], Rarr[GAMG_MAXLEVELS];  PetscReal *coarse_crds = 0, *crds = pc_gamg->m_data;
153*baa4e9faSMark F. Adams   PetscBool isOK;
1545b89ad90SMark F. Adams   for (level=0, Aarr[0] = Pmat; level < GAMG_MAXLEVELS-1; level++ ){
1555b89ad90SMark F. Adams     ierr = MatGetSize( Aarr[level], &M, &N );CHKERRQ(ierr);
156*baa4e9faSMark F. Adams     if( M < npe*10 ) { /* hard wire this for now */
1575b89ad90SMark F. Adams       break;
1585b89ad90SMark F. Adams     }
1595b89ad90SMark F. Adams     level1 = level + 1;
160*baa4e9faSMark F. Adams     PetscPrintf(PETSC_COMM_WORLD,"\t[%d]%s make level %d N=%d\n",0,__FUNCT__,level+1,N);
1610f9369f8SMark F. Adams     ierr = createProlongation( Aarr[level], crds, pc_gamg->m_dim,
162*baa4e9faSMark F. Adams                                &Rarr[level1], &coarse_crds, &isOK );
1635b89ad90SMark F. Adams     CHKERRQ(ierr);
1645b89ad90SMark F. Adams     ierr = PetscFree( crds ); CHKERRQ( ierr );
1655b89ad90SMark F. Adams     crds = coarse_crds;
166*baa4e9faSMark F. Adams     if(level==0) Aarr[0] = Amat; /* use Pmat for finest level setup, but use mat for solver */
167*baa4e9faSMark F. Adams     if( isOK ) {
168*baa4e9faSMark F. Adams       ierr = createCrsOp( Aarr[level], Rarr[level1], &Aarr[level1] ); CHKERRQ(ierr);
169*baa4e9faSMark F. Adams     }
170*baa4e9faSMark F. Adams     else{
171*baa4e9faSMark F. Adams       break;
172*baa4e9faSMark F. Adams     }
1735b89ad90SMark F. Adams   }
1745b89ad90SMark F. Adams   ierr = PetscFree( coarse_crds ); CHKERRQ( ierr );
175*baa4e9faSMark F. Adams PetscPrintf(PETSC_COMM_WORLD,"\t[%d]%s %d levels\n",0,__FUNCT__,level + 1);
1765b89ad90SMark F. Adams   pc_gamg->m_data = 0; /* destroyed coordinate data */
1775b89ad90SMark F. Adams   pc_gamg->m_Nlevels = level + 1;
1785b89ad90SMark F. Adams   fine_level = level;
1795b89ad90SMark F. Adams   ierr = PCMGSetLevels(pc,pc_gamg->m_Nlevels,PETSC_NULL);CHKERRQ(ierr);
1805b89ad90SMark F. Adams 
1815b89ad90SMark F. Adams   /* set default smoothers */
1825b89ad90SMark F. Adams   PetscReal emax = 2.0, emin;
1835b89ad90SMark F. Adams   for (level=1; level<=fine_level; level++){
1845b89ad90SMark F. Adams     KSP smoother; PC subpc;
1855b89ad90SMark F. Adams     ierr = PCMGGetSmoother(pc,level,&smoother);CHKERRQ(ierr);
1865b89ad90SMark F. Adams     ierr = KSPSetType(smoother,KSPCHEBYCHEV);CHKERRQ(ierr);
1870f9369f8SMark F. Adams     emin = emax/10.0; /* fix!!! */
1885b89ad90SMark F. Adams     ierr = KSPGetPC(smoother,&subpc);CHKERRQ(ierr);
1895b89ad90SMark F. Adams     ierr = PCSetType(subpc,PCJACOBI);CHKERRQ(ierr);
1905b89ad90SMark F. Adams     ierr = KSPChebychevSetEigenvalues(smoother, emax, emin);CHKERRQ(ierr); /* need auto !!!!*/
1915b89ad90SMark F. Adams   }
1925b89ad90SMark F. Adams   ierr = PCSetFromOptions_MG(pc); CHKERRQ(ierr); /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */
1935b89ad90SMark F. Adams   {
1945b89ad90SMark F. Adams     PetscBool galerkin;
1955b89ad90SMark F. Adams     ierr = PCMGGetGalerkin( pc,  &galerkin); CHKERRQ(ierr);
1965b89ad90SMark F. Adams     if(galerkin){
1975b89ad90SMark F. Adams       SETERRQ(wcomm,PETSC_ERR_ARG_WRONG, "GAMG does galerkin manually so it must not be used in PC_MG.");
1985b89ad90SMark F. Adams     }
1995b89ad90SMark F. Adams   }
2005b89ad90SMark F. Adams   /* create coarse level and the interpolation between the levels */
2015b89ad90SMark F. Adams   for (level=0,lidx=pc_gamg->m_Nlevels-1; level<fine_level; level++,lidx--){
2025b89ad90SMark F. Adams     level1 = level + 1;
2035b89ad90SMark F. Adams     /* PetscInt MM,NN; */
2045b89ad90SMark F. Adams     /* ierr = MatGetSize( Rarr[lidx], &MM, &NN );CHKERRQ(ierr); */
2055b89ad90SMark F. Adams     /* PetscPrintf(PETSC_COMM_WORLD,"%s Set P(%d,%d) on level %d (%d)\n",__FUNCT__,MM,NN,level1,lidx); */
2065b89ad90SMark F. Adams     ierr = PCMGSetInterpolation(pc,level1,Rarr[lidx]);CHKERRQ(ierr);
2076c237d78SBarry Smith     if(!PETSC_TRUE) {
2085b89ad90SMark F. Adams       PetscViewer        viewer;
2095b89ad90SMark F. Adams       ierr = PetscViewerASCIIOpen(PETSC_COMM_SELF, "Rmat.m", &viewer);  CHKERRQ(ierr);
2105b89ad90SMark F. Adams       ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB);  CHKERRQ(ierr);
211a92563c5SMark F. Adams       ierr = MatView(Rarr[lidx],viewer);CHKERRQ(ierr);
2125b89ad90SMark F. Adams       ierr = PetscViewerDestroy( &viewer );
2135b89ad90SMark F. Adams     }
214a92563c5SMark F. Adams     ierr = MatDestroy( &Rarr[lidx] );  CHKERRQ(ierr);
215a92563c5SMark F. Adams     {
2165b89ad90SMark F. Adams       KSP smoother;
2175b89ad90SMark F. Adams       ierr = PCMGGetSmoother(pc,level,&smoother); CHKERRQ(ierr);
2185b89ad90SMark F. Adams       ierr = KSPSetOperators( smoother, Aarr[lidx], Aarr[lidx], DIFFERENT_NONZERO_PATTERN );
2195b89ad90SMark F. Adams       CHKERRQ(ierr);
220a92563c5SMark F. Adams       ierr = MatDestroy( &Aarr[lidx] );  CHKERRQ(ierr);
221a92563c5SMark F. Adams     }
2225b89ad90SMark F. Adams   }
2235b89ad90SMark F. Adams   { /* fine level (no P) */
2245b89ad90SMark F. Adams     KSP smoother;
2255b89ad90SMark F. Adams     ierr = PCMGGetSmoother(pc,fine_level,&smoother); CHKERRQ(ierr);
2265b89ad90SMark F. Adams     ierr = KSPSetOperators( smoother, Aarr[0], Aarr[0], DIFFERENT_NONZERO_PATTERN );
2275b89ad90SMark F. Adams     CHKERRQ(ierr);
2285b89ad90SMark F. Adams   }
2295b89ad90SMark F. Adams 
2305b89ad90SMark F. Adams   /* setupcalled is set to 0 so that MG is setup from scratch */
2315b89ad90SMark F. Adams   pc->setupcalled = 0;
2325b89ad90SMark F. Adams   ierr = PCSetUp_MG(pc);CHKERRQ(ierr);
2335b89ad90SMark F. Adams   PetscFunctionReturn(0);
2345b89ad90SMark F. Adams }
2355b89ad90SMark F. Adams 
2365b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
2375b89ad90SMark F. Adams /*
2385b89ad90SMark F. Adams    PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner
2395b89ad90SMark F. Adams    that was created with PCCreate_GAMG().
2405b89ad90SMark F. Adams 
2415b89ad90SMark F. Adams    Input Parameter:
2425b89ad90SMark F. Adams .  pc - the preconditioner context
2435b89ad90SMark F. Adams 
2445b89ad90SMark F. Adams    Application Interface Routine: PCDestroy()
2455b89ad90SMark F. Adams */
2465b89ad90SMark F. Adams #undef __FUNCT__
2475b89ad90SMark F. Adams #define __FUNCT__ "PCDestroy_GAMG"
2485b89ad90SMark F. Adams PetscErrorCode PCDestroy_GAMG(PC pc)
2495b89ad90SMark F. Adams {
2505b89ad90SMark F. Adams   PetscErrorCode  ierr;
2515b89ad90SMark F. Adams   PC_MG           *mg = (PC_MG*)pc->data;
2525b89ad90SMark F. Adams   PC_GAMG         *pc_gamg= (PC_GAMG*)mg->innerctx;
2535b89ad90SMark F. Adams 
2545b89ad90SMark F. Adams   PetscFunctionBegin;
2555b89ad90SMark F. Adams   ierr = PCReset_GAMG(pc);CHKERRQ(ierr);
2565b89ad90SMark F. Adams   ierr = PetscFree(pc_gamg);CHKERRQ(ierr);
2575b89ad90SMark F. Adams   ierr = PCDestroy_MG(pc);CHKERRQ(ierr);
2585b89ad90SMark F. Adams   PetscFunctionReturn(0);
2595b89ad90SMark F. Adams }
2605b89ad90SMark F. Adams 
2615b89ad90SMark F. Adams #undef __FUNCT__
2625b89ad90SMark F. Adams #define __FUNCT__ "PCSetFromOptions_GAMG"
2635b89ad90SMark F. Adams PetscErrorCode PCSetFromOptions_GAMG(PC pc)
2645b89ad90SMark F. Adams {
2655b89ad90SMark F. Adams   /* PetscErrorCode  ierr; */
2665b89ad90SMark F. Adams   /* PC_MG           *mg = (PC_MG*)pc->data; */
2675b89ad90SMark F. Adams   /* PC_GAMG         *pc_gamg = (PC_GAMG*)mg->innerctx; */
2685b89ad90SMark F. Adams   /* MPI_Comm        comm = ((PetscObject)pc)->comm; */
2695b89ad90SMark F. Adams 
2705b89ad90SMark F. Adams   PetscFunctionBegin;
2715b89ad90SMark F. Adams   PetscFunctionReturn(0);
2725b89ad90SMark F. Adams }
2735b89ad90SMark F. Adams 
2745b89ad90SMark F. Adams /* -------------------------------------------------------------------------- */
2755b89ad90SMark F. Adams /*
2765b89ad90SMark F. Adams  PCCreate_GAMG - Creates a GAMG preconditioner context, PC_GAMG
2775b89ad90SMark F. Adams 
2785b89ad90SMark F. Adams    Input Parameter:
2795b89ad90SMark F. Adams .  pc - the preconditioner context
2805b89ad90SMark F. Adams 
2815b89ad90SMark F. Adams    Application Interface Routine: PCCreate()
2825b89ad90SMark F. Adams 
2835b89ad90SMark F. Adams   */
2845b89ad90SMark F. Adams  /* MC
2855b89ad90SMark F. Adams      PCGAMG - Use algebraic multigrid preconditioning. This preconditioner requires you provide
2865b89ad90SMark F. Adams        fine grid discretization matrix and coordinates on the fine grid.
2875b89ad90SMark F. Adams 
2885b89ad90SMark F. Adams    Options Database Key:
2895b89ad90SMark F. Adams    Multigrid options(inherited)
2905b89ad90SMark F. Adams +  -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (MGSetCycles)
2915b89ad90SMark F. Adams .  -pc_mg_smoothup <1>: Number of post-smoothing steps (MGSetNumberSmoothUp)
2925b89ad90SMark F. Adams .  -pc_mg_smoothdown <1>: Number of pre-smoothing steps (MGSetNumberSmoothDown)
2935b89ad90SMark F. Adams    -pc_mg_type <multiplicative>: (one of) additive multiplicative full cascade kascade
2945b89ad90SMark F. Adams    GAMG options:
2955b89ad90SMark F. Adams 
2965b89ad90SMark F. Adams    Level: intermediate
2975b89ad90SMark F. Adams   Concepts: multigrid
2985b89ad90SMark F. Adams 
2995b89ad90SMark F. Adams .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType,
3005b89ad90SMark F. Adams            PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), MPSetCycles(), PCMGSetNumberSmoothDown(),
3015b89ad90SMark F. Adams            PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(),
3025b89ad90SMark F. Adams            PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(),
3035b89ad90SMark F. Adams            PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR()
3045b89ad90SMark F. Adams M */
3055b89ad90SMark F. Adams 
3065b89ad90SMark F. Adams EXTERN_C_BEGIN
3075b89ad90SMark F. Adams #undef __FUNCT__
3085b89ad90SMark F. Adams #define __FUNCT__ "PCCreate_GAMG"
3095b89ad90SMark F. Adams PetscErrorCode  PCCreate_GAMG(PC pc)
3105b89ad90SMark F. Adams {
3115b89ad90SMark F. Adams   PetscErrorCode  ierr;
3125b89ad90SMark F. Adams   PC_GAMG         *pc_gamg;
3135b89ad90SMark F. Adams   PC_MG           *mg;
3145b89ad90SMark F. Adams 
3155b89ad90SMark F. Adams   PetscFunctionBegin;
3165b89ad90SMark F. Adams   /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */
3175b89ad90SMark F. Adams   ierr = PCSetType(pc,PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */
3185b89ad90SMark F. Adams   ierr = PetscObjectChangeTypeName((PetscObject)pc,PCGAMG);CHKERRQ(ierr);
3195b89ad90SMark F. Adams 
3205b89ad90SMark F. Adams   /* create a supporting struct and attach it to pc */
3215b89ad90SMark F. Adams   ierr = PetscNewLog(pc,PC_GAMG,&pc_gamg);CHKERRQ(ierr);
3225b89ad90SMark F. Adams   mg = (PC_MG*)pc->data;
3235b89ad90SMark F. Adams   mg->innerctx = pc_gamg;
3245b89ad90SMark F. Adams 
3255b89ad90SMark F. Adams   pc_gamg->m_Nlevels    = -1;
3265b89ad90SMark F. Adams 
3275b89ad90SMark F. Adams   /* overwrite the pointers of PCMG by the functions of PCGAMG */
3285b89ad90SMark F. Adams   pc->ops->setfromoptions = PCSetFromOptions_GAMG;
3295b89ad90SMark F. Adams   pc->ops->setup          = PCSetUp_GAMG;
3305b89ad90SMark F. Adams   pc->ops->reset          = PCReset_GAMG;
3315b89ad90SMark F. Adams   pc->ops->destroy        = PCDestroy_GAMG;
3325b89ad90SMark F. Adams 
3335b89ad90SMark F. Adams   ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc,
3345b89ad90SMark F. Adams 					    "PCSetCoordinates_C",
3355b89ad90SMark F. Adams 					    "PCSetCoordinates_GAMG",
3365b89ad90SMark F. Adams 					    PCSetCoordinates_GAMG);CHKERRQ(ierr);
3375b89ad90SMark F. Adams   PetscFunctionReturn(0);
3385b89ad90SMark F. Adams }
3395b89ad90SMark F. Adams EXTERN_C_END
340