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