1 /* 2 GAMG geometric-algebric multiogrid PC - Mark Adams 2011 3 */ 4 #include <private/pcimpl.h> /*I "petscpc.h" I*/ 5 #include <../src/ksp/pc/impls/mg/mgimpl.h> /*I "petscpcmg.h" I*/ 6 #include <../src/mat/impls/aij/seq/aij.h> 7 #include <../src/mat/impls/aij/mpi/mpiaij.h> 8 #include <assert.h> 9 10 /* Private context for the GAMG preconditioner */ 11 typedef struct gamg_TAG { 12 PetscInt m_dim; 13 PetscInt m_Nlevels; 14 PetscInt m_data_sz; 15 PetscReal *m_data; /* blocked vector of vertex data on fine grid (coordinates) */ 16 } PC_GAMG; 17 18 /* -----------------------------------------------------------------------------*/ 19 #undef __FUNCT__ 20 #define __FUNCT__ "PCReset_GAMG" 21 PetscErrorCode PCReset_GAMG(PC pc) 22 { 23 PetscErrorCode ierr; 24 PC_MG *mg = (PC_MG*)pc->data; 25 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 26 27 PetscFunctionBegin; 28 ierr = PetscFree(pc_gamg->m_data);CHKERRQ(ierr); 29 PetscFunctionReturn(0); 30 } 31 32 /* -------------------------------------------------------------------------- */ 33 /* 34 PCSetCoordinates_GAMG 35 36 Input Parameter: 37 . pc - the preconditioner context 38 */ 39 EXTERN_C_BEGIN 40 #undef __FUNCT__ 41 #define __FUNCT__ "PCSetCoordinates_GAMG" 42 PetscErrorCode PCSetCoordinates_GAMG(PC pc, const int ndm,PetscReal *coords ) 43 { 44 PC_MG *mg = (PC_MG*)pc->data; 45 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 46 PetscErrorCode ierr; 47 PetscInt bs, my0, tt; 48 Mat mat = pc->pmat; 49 PetscInt arrsz; 50 51 PetscFunctionBegin; 52 ierr = MatGetBlockSize( mat, &bs ); CHKERRQ( ierr ); 53 ierr = MatGetOwnershipRange( mat, &my0, &tt ); CHKERRQ(ierr); 54 arrsz = (tt-my0)/bs*ndm; 55 56 // put coordinates 57 if (!pc_gamg->m_data || (pc_gamg->m_data_sz != arrsz)) { 58 ierr = PetscFree( pc_gamg->m_data ); CHKERRQ(ierr); 59 ierr = PetscMalloc(arrsz*sizeof(double), &pc_gamg->m_data ); CHKERRQ(ierr); 60 } 61 62 /* copy data in */ 63 for(tt=0;tt<arrsz;tt++){ 64 pc_gamg->m_data[tt] = coords[tt]; 65 } 66 pc_gamg->m_data_sz = arrsz; 67 pc_gamg->m_dim = ndm; 68 69 PetscFunctionReturn(0); 70 } 71 EXTERN_C_END 72 73 /* -------------------------------------------------------------------------- */ 74 /* 75 partitionLevel 76 77 Input Parameter: 78 . a_Amat_fine - matrix on this fine (k) level 79 . a_dime - 2 or 3 80 In/Output Parameter: 81 . a_P_inout - prolongation operator to the next level (k-1) 82 . a_coarse_crds - coordinates that need to be moved 83 . a_active_proc - number of active procs 84 Output Parameter: 85 . a_Amat_crs - coarse matrix that is created (k-1) 86 */ 87 88 PetscMPIInt n_active_procs; // hack!!!!! 89 90 #undef __FUNCT__ 91 #define __FUNCT__ "partitionLevel" 92 PetscErrorCode partitionLevel( Mat a_Amat_fine, 93 PetscInt a_dim, 94 Mat *a_P_inout, 95 PetscReal **a_coarse_crds, 96 PetscMPIInt *a_active_proc, 97 Mat *a_Amat_crs 98 ) 99 { 100 PetscErrorCode ierr; 101 Mat Amat, Pnew, Pold = *a_P_inout; 102 IS new_indices,isnum; 103 MPI_Comm wcomm = ((PetscObject)a_Amat_fine)->comm; 104 PetscMPIInt nactive,mype,npe; 105 PetscInt Istart,Iend,Istart0,Iend0,ncrs0,ncrs_new,bs=1; /* bs ??? */ 106 107 PetscFunctionBegin; 108 ierr = MPI_Comm_rank(wcomm,&mype);CHKERRQ(ierr); 109 ierr = MPI_Comm_size(wcomm,&npe);CHKERRQ(ierr); 110 /* RAP */ 111 ierr = MatPtAP( a_Amat_fine, Pold, MAT_INITIAL_MATRIX, 2.0, &Amat ); CHKERRQ(ierr); 112 ierr = MatGetOwnershipRange( Amat, &Istart0, &Iend0 ); CHKERRQ(ierr); /* x2 size of 'a_coarse_crds' */ 113 ncrs0 = (Iend0 - Istart0)/bs; 114 115 /* Repartition Amat_{k} and move colums of P^{k}_{k-1} and coordinates accordingly */ 116 { 117 PetscInt neq,N,counts[npe]; 118 IS isnewproc; 119 PetscMPIInt factstart,fact,new_npe,targ_npe; 120 121 ierr = MatGetSize( Amat, &neq, &N );CHKERRQ(ierr); 122 #define MIN_EQ_PROC 100 123 nactive = *a_active_proc; 124 targ_npe = neq/MIN_EQ_PROC; /* hardwire min. number of eq/proc */ 125 if( targ_npe == 0 || neq < 1000 ) new_npe = 1; /* chop coarsest grid */ 126 else if (targ_npe > nactive ) new_npe = nactive; /* no change */ 127 else { 128 new_npe = -9999; 129 factstart = nactive; 130 for(fact=factstart;fact>0;fact--){ /* try to find a better number of procs */ 131 if( nactive%fact==0 && neq/(nactive/fact) > MIN_EQ_PROC ) { 132 new_npe = nactive/fact; 133 } 134 } 135 assert(new_npe != -9999); 136 } 137 *a_active_proc = new_npe; /* output for next time */ 138 139 n_active_procs = new_npe; // hack!!!!! 140 PetscPrintf(PETSC_COMM_WORLD,"\t\t%s n_active_procs=%d\n",__FUNCT__,n_active_procs); 141 { /* partition: get 'isnewproc' */ 142 MatPartitioning mpart; 143 PetscMPIInt factor = npe/new_npe; 144 /* create sub communicator */ 145 MPI_Comm cm,new_comm; 146 int membershipKey = mype % factor; 147 148 ierr = MPI_Comm_split(wcomm, membershipKey, mype, &cm); CHKERRQ(ierr); 149 ierr = PetscCommDuplicate( cm, &new_comm, PETSC_NULL ); CHKERRQ(ierr); 150 ierr = MPI_Comm_free( &cm ); CHKERRQ(ierr); 151 152 /* use 'new_comm' -- todo!!! */ 153 ierr = MatPartitioningCreate( wcomm, &mpart ); CHKERRQ(ierr); 154 ierr = MatPartitioningSetAdjacency( mpart, Amat ); CHKERRQ(ierr); 155 ierr = MatPartitioningSetFromOptions( mpart ); CHKERRQ(ierr); 156 ierr = MatPartitioningApply( mpart, &isnewproc ); CHKERRQ(ierr); 157 ierr = MatPartitioningDestroy( &mpart ); CHKERRQ(ierr); 158 159 ierr = MPI_Comm_free( &new_comm ); CHKERRQ(ierr); 160 161 /* spread partitioning across machine - probably the right thing to do but machine specific */ 162 /* how to do this??? scale IS by 'factor' */ 163 } 164 165 /* 166 Create an index set from the isnewproc index set to indicate the mapping TO 167 */ 168 ierr = ISPartitioningToNumbering( isnewproc, &isnum ); CHKERRQ(ierr); 169 /* 170 Determine how many elements are assigned to each processor 171 */ 172 ierr = ISPartitioningCount( isnewproc, npe, counts ); CHKERRQ(ierr); 173 ierr = ISDestroy( &isnewproc ); CHKERRQ(ierr); 174 ncrs_new = counts[mype]; 175 PetscPrintf(PETSC_COMM_SELF,"\t[%d]%s local equations = %d\n",mype,__FUNCT__,ncrs_new); 176 } 177 { /* Create a vector to contain the newly ordered element information */ 178 const PetscInt *idx; 179 PetscInt i,j,k; 180 IS isscat; 181 PetscScalar *array; 182 Vec src_crd, dest_crd; 183 PetscReal *coords = *a_coarse_crds; 184 VecScatter vecscat; 185 186 ierr = VecCreate( wcomm, &dest_crd ); 187 ierr = VecSetSizes( dest_crd, a_dim*ncrs_new, PETSC_DECIDE ); CHKERRQ(ierr); 188 ierr = VecSetFromOptions( dest_crd ); CHKERRQ(ierr); /*funny vector-get global options?*/ 189 /* 190 There are three data items per cell (element), the integer vertex numbers of its three 191 coordinates (we convert to double to use the scatter) (one can think 192 of the vectors of having a block size of 3, then there is one index in idx[] for each element) 193 */ 194 { 195 PetscInt tidx[ncrs0*a_dim]; 196 ierr = ISGetIndices( isnum, &idx ); CHKERRQ(ierr); 197 for (i=0,j=0; i<ncrs0; i++) { 198 for (k=0; k<a_dim; k++) tidx[j++] = idx[i]*a_dim + k; 199 } 200 ierr = ISCreateGeneral( wcomm, a_dim*ncrs0, tidx, PETSC_COPY_VALUES, &isscat ); 201 CHKERRQ(ierr); 202 ierr = ISRestoreIndices( isnum, &idx ); CHKERRQ(ierr); 203 } 204 /* 205 Create a vector to contain the original vertex information for each element 206 */ 207 ierr = VecCreateSeq( PETSC_COMM_SELF, a_dim*ncrs0, &src_crd ); CHKERRQ(ierr); 208 ierr = VecGetArray( src_crd, &array ); CHKERRQ(ierr); 209 for (i=0; i<a_dim*ncrs0; i++) array[i] = coords[i]; 210 ierr = VecRestoreArray( src_crd, &array ); CHKERRQ(ierr); 211 /* 212 Scatter the element vertex information (still in the original vertex ordering) 213 to the correct processor 214 */ 215 ierr = VecScatterCreate( src_crd, PETSC_NULL, dest_crd, isscat, &vecscat); 216 CHKERRQ(ierr); 217 ierr = ISDestroy( &isscat ); CHKERRQ(ierr); 218 ierr = VecScatterBegin(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 219 ierr = VecScatterEnd(vecscat,src_crd,dest_crd,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 220 ierr = VecScatterDestroy( &vecscat ); CHKERRQ(ierr); 221 ierr = VecDestroy( &src_crd ); CHKERRQ(ierr); 222 /* 223 Put the element vertex data into a new allocation of the gdata->ele 224 */ 225 ierr = PetscFree( *a_coarse_crds ); CHKERRQ(ierr); 226 ierr = PetscMalloc( a_dim*ncrs_new*sizeof(PetscReal), a_coarse_crds ); CHKERRQ(ierr); 227 coords = *a_coarse_crds; /* convient */ 228 ierr = VecGetArray( dest_crd, &array ); CHKERRQ(ierr); 229 for (i=0; i<a_dim*ncrs_new; i++) coords[i] = PetscRealPart(array[i]); 230 ierr = VecRestoreArray( dest_crd, &array ); CHKERRQ(ierr); 231 ierr = VecDestroy( &dest_crd ); CHKERRQ(ierr); 232 } 233 /* 234 Invert for MatGetSubMatrix 235 */ 236 ierr = ISInvertPermutation( isnum, ncrs_new, &new_indices ); CHKERRQ(ierr); 237 ierr = ISSort( new_indices ); CHKERRQ(ierr); /* is this needed? */ 238 ierr = ISDestroy( &isnum ); CHKERRQ(ierr); 239 /* A_crs output */ 240 ierr = MatGetSubMatrix( Amat, new_indices, new_indices, MAT_INITIAL_MATRIX, a_Amat_crs ); 241 CHKERRQ(ierr); 242 ierr = MatDestroy( &Amat ); CHKERRQ(ierr); 243 Amat = *a_Amat_crs; 244 /* prolongator */ 245 ierr = MatGetOwnershipRange( Pold, &Istart, &Iend ); CHKERRQ(ierr); 246 { 247 IS findices; 248 ierr = ISCreateStride(wcomm,Iend-Istart,Istart,1,&findices); CHKERRQ(ierr); 249 ierr = MatGetSubMatrix( Pold, findices, new_indices, MAT_INITIAL_MATRIX, &Pnew ); 250 CHKERRQ(ierr); 251 ierr = ISDestroy( &findices ); CHKERRQ(ierr); 252 } 253 ierr = MatDestroy( a_P_inout ); CHKERRQ(ierr); 254 *a_P_inout = Pnew; /* output */ 255 256 ierr = ISDestroy( &new_indices ); CHKERRQ(ierr); 257 258 PetscFunctionReturn(0); 259 } 260 261 /* -------------------------------------------------------------------------- */ 262 /* 263 PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner 264 by setting data structures and options. 265 266 Input Parameter: 267 . pc - the preconditioner context 268 269 Application Interface Routine: PCSetUp() 270 271 Notes: 272 The interface routine PCSetUp() is not usually called directly by 273 the user, but instead is called by PCApply() if necessary. 274 */ 275 extern PetscErrorCode PCSetFromOptions_MG(PC); 276 extern PetscErrorCode PCReset_MG(PC); 277 extern PetscErrorCode createProlongation( Mat, PetscReal [], const PetscInt, 278 Mat *, PetscReal **, PetscBool *a_isOK ); 279 #undef __FUNCT__ 280 #define __FUNCT__ "PCSetUp_GAMG" 281 PetscErrorCode PCSetUp_GAMG( PC pc ) 282 { 283 PetscErrorCode ierr; 284 PC_MG *mg = (PC_MG*)pc->data; 285 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 286 Mat Amat = pc->mat, Pmat = pc->pmat; 287 PetscBool isSeq, isMPI; 288 PetscInt fine_level, level, level1, M, N, bs, lidx; 289 MPI_Comm wcomm = ((PetscObject)pc)->comm; 290 PetscMPIInt mype,npe,nactivepe; 291 292 PetscFunctionBegin; 293 ierr = MPI_Comm_rank(wcomm,&mype);CHKERRQ(ierr); 294 ierr = MPI_Comm_size(wcomm,&npe);CHKERRQ(ierr); 295 if (pc->setupcalled){ 296 /* no state data in GAMG to destroy (now) */ 297 ierr = PCReset_MG(pc); CHKERRQ(ierr); 298 } 299 if (!pc_gamg->m_data) SETERRQ(wcomm,PETSC_ERR_SUP,"PCSetUp_GAMG called before PCSetCoordinates"); 300 /* setup special features of PCGAMG */ 301 ierr = PetscTypeCompare((PetscObject)Amat, MATSEQAIJ, &isSeq);CHKERRQ(ierr); 302 ierr = PetscTypeCompare((PetscObject)Amat, MATMPIAIJ, &isMPI);CHKERRQ(ierr); 303 if (isMPI) { 304 } else if (isSeq) { 305 } 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); 306 307 /* GAMG requires input of fine-grid matrix. It determines nlevels. */ 308 ierr = MatGetBlockSize( Amat, &bs ); CHKERRQ(ierr); 309 if(bs!=1) SETERRQ1(wcomm,PETSC_ERR_ARG_WRONG, "GAMG only supports scalar prblems bs = '%d'.",bs); 310 311 /* Get A_i and R_i */ 312 ierr = MatGetSize( Amat, &M, &N );CHKERRQ(ierr); 313 PetscPrintf(PETSC_COMM_WORLD,"[%d]%s level %d N=%d\n",0,__FUNCT__,0,N); 314 #define GAMG_MAXLEVELS 20 315 Mat Aarr[GAMG_MAXLEVELS], Rarr[GAMG_MAXLEVELS]; PetscReal *coarse_crds = 0, *crds = pc_gamg->m_data; 316 PetscBool isOK; 317 for (level=0, Aarr[0] = Pmat, nactivepe = npe; level < GAMG_MAXLEVELS-1; level++ ){ 318 if( nactivepe == 1 ) { 319 break; 320 } 321 level1 = level + 1; 322 ierr = createProlongation( Aarr[level], crds, pc_gamg->m_dim, 323 &Rarr[level1], &coarse_crds, &isOK ); 324 CHKERRQ(ierr); 325 ierr = PetscFree( crds ); CHKERRQ( ierr ); 326 if(level==0) Aarr[0] = Amat; /* use Pmat for finest level setup, but use mat for solver */ 327 if( isOK ) { 328 ierr = partitionLevel( Aarr[level], pc_gamg->m_dim, &Rarr[level1], &coarse_crds, &nactivepe, &Aarr[level1] ); 329 CHKERRQ(ierr); 330 ierr = MatGetSize( Aarr[level1], &M, &N );CHKERRQ(ierr); 331 PetscPrintf(PETSC_COMM_WORLD,"[%d]%s done level %d N=%d, active pe = %d\n",0,__FUNCT__,level1,N,nactivepe); 332 } 333 else{ 334 break; 335 } 336 crds = coarse_crds; 337 } 338 ierr = PetscFree( coarse_crds ); CHKERRQ( ierr ); 339 PetscPrintf(PETSC_COMM_WORLD,"\t[%d]%s %d levels\n",0,__FUNCT__,level + 1); 340 pc_gamg->m_data = 0; /* destroyed coordinate data */ 341 pc_gamg->m_Nlevels = level + 1; 342 fine_level = level; 343 ierr = PCMGSetLevels(pc,pc_gamg->m_Nlevels,PETSC_NULL);CHKERRQ(ierr); 344 345 /* set default smoothers */ 346 PetscReal emax = 2.0, emin; 347 for (level=1; level<=fine_level; level++){ 348 KSP smoother; PC subpc; 349 ierr = PCMGGetSmoother(pc,level,&smoother);CHKERRQ(ierr); 350 ierr = KSPSetType(smoother,KSPCHEBYCHEV);CHKERRQ(ierr); 351 emin = emax/10.0; /* fix!!! */ 352 ierr = KSPGetPC(smoother,&subpc);CHKERRQ(ierr); 353 ierr = PCSetType(subpc,PCJACOBI);CHKERRQ(ierr); 354 ierr = KSPChebychevSetEigenvalues(smoother, emax, emin);CHKERRQ(ierr); /* need auto !!!!*/ 355 } 356 ierr = PCSetFromOptions_MG(pc); CHKERRQ(ierr); /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */ 357 { 358 PetscBool galerkin; 359 ierr = PCMGGetGalerkin( pc, &galerkin); CHKERRQ(ierr); 360 if(galerkin){ 361 SETERRQ(wcomm,PETSC_ERR_ARG_WRONG, "GAMG does galerkin manually so it must not be used in PC_MG."); 362 } 363 } 364 /* create coarse level and the interpolation between the levels */ 365 for (level=0,lidx=pc_gamg->m_Nlevels-1; level<fine_level; level++,lidx--){ 366 level1 = level + 1; 367 ierr = PCMGSetInterpolation(pc,level1,Rarr[lidx]);CHKERRQ(ierr); 368 if(!PETSC_TRUE) { 369 PetscViewer viewer; char fname[32]; 370 sprintf(fname,"Amat_%d.m",lidx); 371 ierr = PetscViewerASCIIOpen( wcomm, fname, &viewer ); CHKERRQ(ierr); 372 ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB); CHKERRQ(ierr); 373 ierr = MatView( Aarr[lidx], viewer ); CHKERRQ(ierr); 374 ierr = PetscViewerDestroy( &viewer ); 375 } 376 ierr = MatDestroy( &Rarr[lidx] ); CHKERRQ(ierr); 377 { 378 KSP smoother; 379 ierr = PCMGGetSmoother(pc,level,&smoother); CHKERRQ(ierr); 380 ierr = KSPSetOperators( smoother, Aarr[lidx], Aarr[lidx], DIFFERENT_NONZERO_PATTERN ); 381 CHKERRQ(ierr); 382 ierr = MatDestroy( &Aarr[lidx] ); CHKERRQ(ierr); 383 } 384 } 385 { /* fine level (no P) */ 386 KSP smoother; 387 ierr = PCMGGetSmoother(pc,fine_level,&smoother); CHKERRQ(ierr); 388 ierr = KSPSetOperators( smoother, Aarr[0], Aarr[0], DIFFERENT_NONZERO_PATTERN ); 389 CHKERRQ(ierr); 390 } 391 392 /* setupcalled is set to 0 so that MG is setup from scratch */ 393 pc->setupcalled = 0; 394 ierr = PCSetUp_MG(pc);CHKERRQ(ierr); 395 PetscFunctionReturn(0); 396 } 397 398 /* -------------------------------------------------------------------------- */ 399 /* 400 PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner 401 that was created with PCCreate_GAMG(). 402 403 Input Parameter: 404 . pc - the preconditioner context 405 406 Application Interface Routine: PCDestroy() 407 */ 408 #undef __FUNCT__ 409 #define __FUNCT__ "PCDestroy_GAMG" 410 PetscErrorCode PCDestroy_GAMG(PC pc) 411 { 412 PetscErrorCode ierr; 413 PC_MG *mg = (PC_MG*)pc->data; 414 PC_GAMG *pc_gamg= (PC_GAMG*)mg->innerctx; 415 416 PetscFunctionBegin; 417 ierr = PCReset_GAMG(pc);CHKERRQ(ierr); 418 ierr = PetscFree(pc_gamg);CHKERRQ(ierr); 419 ierr = PCDestroy_MG(pc);CHKERRQ(ierr); 420 PetscFunctionReturn(0); 421 } 422 423 #undef __FUNCT__ 424 #define __FUNCT__ "PCSetFromOptions_GAMG" 425 PetscErrorCode PCSetFromOptions_GAMG(PC pc) 426 { 427 /* PetscErrorCode ierr; */ 428 /* PC_MG *mg = (PC_MG*)pc->data; */ 429 /* PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; */ 430 /* MPI_Comm comm = ((PetscObject)pc)->comm; */ 431 432 PetscFunctionBegin; 433 PetscFunctionReturn(0); 434 } 435 436 /* -------------------------------------------------------------------------- */ 437 /* 438 PCCreate_GAMG - Creates a GAMG preconditioner context, PC_GAMG 439 440 Input Parameter: 441 . pc - the preconditioner context 442 443 Application Interface Routine: PCCreate() 444 445 */ 446 /* MC 447 PCGAMG - Use algebraic multigrid preconditioning. This preconditioner requires you provide 448 fine grid discretization matrix and coordinates on the fine grid. 449 450 Options Database Key: 451 Multigrid options(inherited) 452 + -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (MGSetCycles) 453 . -pc_mg_smoothup <1>: Number of post-smoothing steps (MGSetNumberSmoothUp) 454 . -pc_mg_smoothdown <1>: Number of pre-smoothing steps (MGSetNumberSmoothDown) 455 -pc_mg_type <multiplicative>: (one of) additive multiplicative full cascade kascade 456 GAMG options: 457 458 Level: intermediate 459 Concepts: multigrid 460 461 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, 462 PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), MPSetCycles(), PCMGSetNumberSmoothDown(), 463 PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(), 464 PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(), 465 PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR() 466 M */ 467 468 EXTERN_C_BEGIN 469 #undef __FUNCT__ 470 #define __FUNCT__ "PCCreate_GAMG" 471 PetscErrorCode PCCreate_GAMG(PC pc) 472 { 473 PetscErrorCode ierr; 474 PC_GAMG *pc_gamg; 475 PC_MG *mg; 476 477 PetscFunctionBegin; 478 /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */ 479 ierr = PCSetType(pc,PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */ 480 ierr = PetscObjectChangeTypeName((PetscObject)pc,PCGAMG);CHKERRQ(ierr); 481 482 /* create a supporting struct and attach it to pc */ 483 ierr = PetscNewLog(pc,PC_GAMG,&pc_gamg);CHKERRQ(ierr); 484 mg = (PC_MG*)pc->data; 485 mg->innerctx = pc_gamg; 486 487 pc_gamg->m_Nlevels = -1; 488 489 /* overwrite the pointers of PCMG by the functions of PCGAMG */ 490 pc->ops->setfromoptions = PCSetFromOptions_GAMG; 491 pc->ops->setup = PCSetUp_GAMG; 492 pc->ops->reset = PCReset_GAMG; 493 pc->ops->destroy = PCDestroy_GAMG; 494 495 ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc, 496 "PCSetCoordinates_C", 497 "PCSetCoordinates_GAMG", 498 PCSetCoordinates_GAMG);CHKERRQ(ierr); 499 PetscFunctionReturn(0); 500 } 501 EXTERN_C_END 502