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 9 /* Private context for the GAMG preconditioner */ 10 typedef struct gamg_TAG { 11 gamg_TAG() : m_dim(0), m_Nlevels(-1), m_data(0) 12 {} 13 PetscInt m_dim; 14 PetscInt m_Nlevels; 15 PetscInt m_data_sz; 16 PetscReal *m_data; /* blocked vector of vertex data on fine grid (coordinates) */ 17 } PC_GAMG; 18 19 /* -----------------------------------------------------------------------------*/ 20 #undef __FUNCT__ 21 #define __FUNCT__ "PCReset_GAMG" 22 PetscErrorCode PCReset_GAMG(PC pc) 23 { 24 PetscErrorCode ierr; 25 PC_MG *mg = (PC_MG*)pc->data; 26 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 27 28 PetscFunctionBegin; 29 if (pc_gamg->m_data) { 30 ierr = PetscFree(pc_gamg->m_data);CHKERRQ(ierr); 31 pc_gamg->m_data = 0; 32 } 33 PetscFunctionReturn(0); 34 } 35 36 /* -------------------------------------------------------------------------- */ 37 /* 38 PCSetCoordinates_GAMG 39 40 Input Parameter: 41 . pc - the preconditioner context 42 */ 43 #undef __FUNCT__ 44 #define __FUNCT__ "PCSetCoordinates_GAMG" 45 PetscErrorCode PCSetCoordinates_GAMG( PC pc, const int ndm, 46 PetscReal *coords ) 47 { 48 PC_MG *mg = (PC_MG*)pc->data; 49 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 50 PetscErrorCode ierr; PetscInt bs, my0, tt; 51 52 PetscFunctionBegin; 53 Mat mat = pc->pmat; 54 ierr = MatGetBlockSize( mat, &bs ); CHKERRQ( ierr ); 55 ierr = MatGetOwnershipRange( mat, &my0, &tt ); CHKERRQ(ierr); 56 const PetscInt arrsz = (tt-my0)/bs*ndm; 57 58 // put coordinates 59 if ( pc_gamg->m_data==0 || pc_gamg->m_data_sz != arrsz ) { 60 if( pc_gamg->m_data != 0 ) { 61 ierr = PetscFree( pc_gamg->m_data ); CHKERRQ(ierr); 62 } 63 ierr = PetscMalloc( (arrsz+1)*sizeof(double), &pc_gamg->m_data ); CHKERRQ(ierr); 64 } 65 66 /* copy data in */ 67 for(tt=0;tt<arrsz;tt++){ 68 pc_gamg->m_data[tt] = coords[tt]; 69 } 70 pc_gamg->m_data_sz = arrsz; 71 pc_gamg->m_dim = ndm; 72 73 PetscFunctionReturn(0); 74 } 75 76 /* -------------------------------------------------------------------------- */ 77 /* 78 createCrsOp 79 80 Input Parameter: 81 . Amat - matrix on this fine level 82 . P_out - prolongation operator to the next level 83 . Acrs - coarse matrix that is created 84 */ 85 #undef __FUNCT__ 86 #define __FUNCT__ "createCrsOp" 87 PetscErrorCode createCrsOp( Mat Amat, Mat P_inout, Mat *Acrs ) 88 { 89 PetscErrorCode ierr; 90 91 PetscFunctionBegin; 92 Mat H; 93 ierr = MatPtAP( Amat, P_inout, MAT_INITIAL_MATRIX, 2.0, &H); CHKERRQ(ierr); 94 95 /* need to repartition H and move colums of P accordingly */ 96 97 98 *Acrs = H; 99 100 PetscFunctionReturn(0); 101 } 102 103 /* -------------------------------------------------------------------------- */ 104 /* 105 PCSetUp_GAMG - Prepares for the use of the GAMG preconditioner 106 by setting data structures and options. 107 108 Input Parameter: 109 . pc - the preconditioner context 110 111 Application Interface Routine: PCSetUp() 112 113 Notes: 114 The interface routine PCSetUp() is not usually called directly by 115 the user, but instead is called by PCApply() if necessary. 116 */ 117 extern PetscErrorCode PCSetFromOptions_MG(PC); 118 extern PetscErrorCode PCReset_MG(PC); 119 extern PetscErrorCode createProlongation( Mat, Mat *, PetscReal [], PetscReal **, const PetscInt); 120 #undef __FUNCT__ 121 #define __FUNCT__ "PCSetUp_GAMG" 122 PetscErrorCode PCSetUp_GAMG( PC pc ) 123 { 124 PetscErrorCode ierr; 125 PC_MG *mg = (PC_MG*)pc->data; 126 PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; 127 Mat Amat = pc->mat, Pmat = pc->pmat; 128 PetscBool isSeq, isMPI; 129 PetscInt fine_level, level, level1, M, N, bs, lidx; 130 MPI_Comm wcomm = ((PetscObject)pc)->comm; 131 132 PetscFunctionBegin; 133 if (pc->setupcalled){ 134 /* no state data in GAMG to destroy (now) */ 135 ierr = PCReset_MG(pc);CHKERRQ(ierr); 136 } 137 if ( pc_gamg->m_data==0 ) { 138 SETERRQ(wcomm,PETSC_ERR_SUP,"PCSetUp_GAMG called before PCSetCoordinates"); 139 } 140 /* setup special features of PCGAMG */ 141 ierr = PetscTypeCompare((PetscObject) Amat, MATSEQAIJ, &isSeq);CHKERRQ(ierr); 142 ierr = PetscTypeCompare((PetscObject) Amat, MATMPIAIJ, &isMPI);CHKERRQ(ierr); 143 if (isMPI) { 144 } else if (isSeq) { 145 } 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); 146 147 /* GAMG requires input of fine-grid matrix. It determines nlevels. */ 148 ierr = MatGetSize( Amat, &M, &N );CHKERRQ(ierr); 149 ierr = MatGetBlockSize( Amat, &bs ); CHKERRQ(ierr); 150 if(bs!=1) SETERRQ1(wcomm,PETSC_ERR_ARG_WRONG, "GAMG only supports scalar prblems bs = '%d'.",bs); 151 152 /* Get A_i and R_i */ 153 #define GAMG_MAXLEVELS 10 154 Mat Aarr[GAMG_MAXLEVELS], Rarr[GAMG_MAXLEVELS]; PetscReal *coarse_crds = 0, *crds = pc_gamg->m_data; 155 for (level=0, Aarr[0] = Pmat; level < GAMG_MAXLEVELS-1; level++ ){ 156 ierr = MatGetSize( Aarr[level], &M, &N );CHKERRQ(ierr); 157 if( M < 100 ) { /* hard wire this for now */ 158 break; 159 } 160 level1 = level + 1; 161 ierr = createProlongation( Aarr[level], &Rarr[level1], crds, &coarse_crds, pc_gamg->m_dim ); 162 CHKERRQ(ierr); 163 if(level==0) Aarr[0] = Amat; /* use Pmat for finest level setup, but use mat for solver */ 164 ierr = createCrsOp( Aarr[level], Rarr[level1], &Aarr[level1] ); CHKERRQ(ierr); 165 ierr = PetscFree( crds ); CHKERRQ( ierr ); 166 crds = coarse_crds; 167 } 168 if( coarse_crds != 0 ) { 169 ierr = PetscFree( coarse_crds ); CHKERRQ( ierr ); 170 } 171 pc_gamg->m_data = 0; /* destroyed coordinate data */ 172 pc_gamg->m_Nlevels = level + 1; 173 fine_level = level; 174 ierr = PCMGSetLevels(pc,pc_gamg->m_Nlevels,PETSC_NULL);CHKERRQ(ierr); 175 176 /* set default smoothers */ 177 PetscReal emax = 2.0, emin; 178 for (level=1; level<=fine_level; level++){ 179 KSP smoother; PC subpc; 180 ierr = PCMGGetSmoother(pc,level,&smoother);CHKERRQ(ierr); 181 ierr = KSPSetType(smoother,KSPCHEBYCHEV);CHKERRQ(ierr); 182 emin = emax/5.0; 183 ierr = KSPGetPC(smoother,&subpc);CHKERRQ(ierr); 184 ierr = PCSetType(subpc,PCJACOBI);CHKERRQ(ierr); 185 ierr = KSPChebychevSetEigenvalues(smoother, emax, emin);CHKERRQ(ierr); /* need auto !!!!*/ 186 } 187 ierr = PCSetFromOptions_MG(pc); CHKERRQ(ierr); /* should be called in PCSetFromOptions_GAMG(), but cannot be called prior to PCMGSetLevels() */ 188 { 189 PetscBool galerkin; 190 ierr = PCMGGetGalerkin( pc, &galerkin); CHKERRQ(ierr); 191 if(galerkin){ 192 SETERRQ(wcomm,PETSC_ERR_ARG_WRONG, "GAMG does galerkin manually so it must not be used in PC_MG."); 193 } 194 } 195 /* create coarse level and the interpolation between the levels */ 196 for (level=0,lidx=pc_gamg->m_Nlevels-1; level<fine_level; level++,lidx--){ 197 level1 = level + 1; 198 /* PetscInt MM,NN; */ 199 /* ierr = MatGetSize( Rarr[lidx], &MM, &NN );CHKERRQ(ierr); */ 200 /* PetscPrintf(PETSC_COMM_WORLD,"%s Set P(%d,%d) on level %d (%d)\n",__FUNCT__,MM,NN,level1,lidx); */ 201 ierr = PCMGSetInterpolation(pc,level1,Rarr[lidx]);CHKERRQ(ierr); 202 if(!true) { 203 PetscViewer viewer; 204 ierr = PetscViewerASCIIOpen(PETSC_COMM_SELF, "Rmat.m", &viewer); CHKERRQ(ierr); 205 ierr = PetscViewerSetFormat( viewer, PETSC_VIEWER_ASCII_MATLAB); CHKERRQ(ierr); 206 ierr = MatView(Rarr[level1],viewer);CHKERRQ(ierr); 207 ierr = PetscViewerDestroy( &viewer ); 208 } 209 KSP smoother; 210 ierr = PCMGGetSmoother(pc,level,&smoother); CHKERRQ(ierr); 211 ierr = KSPSetOperators( smoother, Aarr[lidx], Aarr[lidx], DIFFERENT_NONZERO_PATTERN ); 212 CHKERRQ(ierr); 213 /* ierr = MatGetSize( Aarr[lidx], &MM, &NN ); CHKERRQ(ierr); */ 214 /* PetscPrintf(PETSC_COMM_WORLD,"%s Set A(%d,%d) on level %d (%d)\n",__FUNCT__,MM,NN,level,lidx); */ 215 } 216 { /* fine level (no P) */ 217 KSP smoother; 218 ierr = PCMGGetSmoother(pc,fine_level,&smoother); CHKERRQ(ierr); 219 ierr = KSPSetOperators( smoother, Aarr[0], Aarr[0], DIFFERENT_NONZERO_PATTERN ); 220 CHKERRQ(ierr); 221 /* PetscInt MM,NN; */ 222 /* ierr = MatGetSize( Aarr[0], &MM, &NN );CHKERRQ(ierr); */ 223 /* PetscPrintf(PETSC_COMM_WORLD,"%s Set A(%d,%d) on level %d (%d)\n",__FUNCT__,MM,NN,fine_level,0); */ 224 } 225 226 /* setupcalled is set to 0 so that MG is setup from scratch */ 227 pc->setupcalled = 0; 228 ierr = PCSetUp_MG(pc);CHKERRQ(ierr); 229 PetscFunctionReturn(0); 230 } 231 232 /* -------------------------------------------------------------------------- */ 233 /* 234 PCDestroy_GAMG - Destroys the private context for the GAMG preconditioner 235 that was created with PCCreate_GAMG(). 236 237 Input Parameter: 238 . pc - the preconditioner context 239 240 Application Interface Routine: PCDestroy() 241 */ 242 #undef __FUNCT__ 243 #define __FUNCT__ "PCDestroy_GAMG" 244 PetscErrorCode PCDestroy_GAMG(PC pc) 245 { 246 PetscErrorCode ierr; 247 PC_MG *mg = (PC_MG*)pc->data; 248 PC_GAMG *pc_gamg= (PC_GAMG*)mg->innerctx; 249 250 PetscFunctionBegin; 251 ierr = PCReset_GAMG(pc);CHKERRQ(ierr); 252 if (pc_gamg->m_data) { 253 ierr = PetscFree(pc_gamg->m_data);CHKERRQ(ierr); 254 } 255 ierr = PetscFree(pc_gamg);CHKERRQ(ierr); 256 ierr = PCDestroy_MG(pc);CHKERRQ(ierr); 257 PetscFunctionReturn(0); 258 } 259 260 #undef __FUNCT__ 261 #define __FUNCT__ "PCSetFromOptions_GAMG" 262 PetscErrorCode PCSetFromOptions_GAMG(PC pc) 263 { 264 /* PetscErrorCode ierr; */ 265 /* PC_MG *mg = (PC_MG*)pc->data; */ 266 /* PC_GAMG *pc_gamg = (PC_GAMG*)mg->innerctx; */ 267 /* MPI_Comm comm = ((PetscObject)pc)->comm; */ 268 269 PetscFunctionBegin; 270 PetscFunctionReturn(0); 271 } 272 273 /* -------------------------------------------------------------------------- */ 274 /* 275 PCCreate_GAMG - Creates a GAMG preconditioner context, PC_GAMG 276 277 Input Parameter: 278 . pc - the preconditioner context 279 280 Application Interface Routine: PCCreate() 281 282 */ 283 /* MC 284 PCGAMG - Use algebraic multigrid preconditioning. This preconditioner requires you provide 285 fine grid discretization matrix and coordinates on the fine grid. 286 287 Options Database Key: 288 Multigrid options(inherited) 289 + -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (MGSetCycles) 290 . -pc_mg_smoothup <1>: Number of post-smoothing steps (MGSetNumberSmoothUp) 291 . -pc_mg_smoothdown <1>: Number of pre-smoothing steps (MGSetNumberSmoothDown) 292 -pc_mg_type <multiplicative>: (one of) additive multiplicative full cascade kascade 293 GAMG options: 294 295 Level: intermediate 296 Concepts: multigrid 297 298 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, 299 PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), MPSetCycles(), PCMGSetNumberSmoothDown(), 300 PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(), 301 PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(), 302 PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR() 303 M */ 304 305 EXTERN_C_BEGIN 306 #undef __FUNCT__ 307 #define __FUNCT__ "PCCreate_GAMG" 308 PetscErrorCode PCCreate_GAMG(PC pc) 309 { 310 PetscErrorCode ierr; 311 PC_GAMG *pc_gamg; 312 PC_MG *mg; 313 314 PetscFunctionBegin; 315 /* PCGAMG is an inherited class of PCMG. Initialize pc as PCMG */ 316 ierr = PCSetType(pc,PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */ 317 ierr = PetscObjectChangeTypeName((PetscObject)pc,PCGAMG);CHKERRQ(ierr); 318 319 /* create a supporting struct and attach it to pc */ 320 ierr = PetscNewLog(pc,PC_GAMG,&pc_gamg);CHKERRQ(ierr); 321 mg = (PC_MG*)pc->data; 322 mg->innerctx = pc_gamg; 323 324 pc_gamg->m_data = 0; 325 pc_gamg->m_Nlevels = -1; 326 327 /* overwrite the pointers of PCMG by the functions of PCGAMG */ 328 pc->ops->setfromoptions = PCSetFromOptions_GAMG; 329 pc->ops->setup = PCSetUp_GAMG; 330 pc->ops->reset = PCReset_GAMG; 331 pc->ops->destroy = PCDestroy_GAMG; 332 333 ierr = PetscObjectComposeFunctionDynamic( (PetscObject)pc, 334 "PCSetCoordinates_C", 335 "PCSetCoordinates_GAMG", 336 PCSetCoordinates_GAMG);CHKERRQ(ierr); 337 PetscFunctionReturn(0); 338 } 339 EXTERN_C_END 340