1 2 /* 3 Provides an interface to the ML 3.0 smoothed Aggregation 4 */ 5 #include "src/ksp/pc/pcimpl.h" /*I "petscpc.h" I*/ 6 #include "src/ksp/pc/impls/mg/mgimpl.h" /*I "petscmg.h" I*/ 7 #include "src/mat/impls/aij/seq/aij.h" 8 #include "src/mat/impls/aij/mpi/mpiaij.h" 9 EXTERN_C_BEGIN 10 #include <math.h> 11 #include "ml_include.h" 12 EXTERN_C_END 13 14 /* The context (data structure) at each grid level */ 15 typedef struct { 16 Vec x,b,r; /* global vectors */ 17 Mat A,P,R; 18 KSP ksp; 19 } GridCtx; 20 21 /* The context used to input PETSc matrix into ML at fine grid */ 22 typedef struct { 23 Mat A,Aloc; 24 Vec x,y; 25 ML_Operator *mlmat; 26 PetscScalar *pwork; /* tmp array used by PetscML_comm() */ 27 } FineGridCtx; 28 29 /* The context associates a ML matrix with a PETSc shell matrix */ 30 typedef struct { 31 Mat A; /* PETSc shell matrix associated with mlmat */ 32 ML_Operator *mlmat; /* ML matrix assorciated with A */ 33 Vec y; 34 } Mat_MLShell; 35 36 /* Private context for the ML preconditioner */ 37 typedef struct { 38 ML *ml_object; 39 ML_Aggregate *agg_object; 40 GridCtx *gridctx; 41 FineGridCtx *PetscMLdata; 42 PetscInt fine_level,MaxNlevels,MaxCoarseSize,CoarsenScheme; 43 PetscReal Threshold,DampingFactor; 44 PetscTruth SpectralNormScheme_Anorm; 45 PetscMPIInt size; 46 PetscErrorCode (*PCSetUp)(PC); 47 PetscErrorCode (*PCDestroy)(PC); 48 } PC_ML; 49 extern int PetscML_getrow(void *ML_data,int N_requested_rows,int requested_rows[], 50 int allocated_space,int columns[],double values[],int row_lengths[]); 51 extern int PetscML_matvec(void *ML_data, int in_length, double p[], int out_length,double ap[]); 52 extern int PetscML_comm(double x[], void *ML_data); 53 extern PetscErrorCode MatMult_ML(Mat,Vec,Vec); 54 extern PetscErrorCode MatMultAdd_ML(Mat,Vec,Vec,Vec); 55 extern PetscErrorCode MatConvert_MPIAIJ_ML(Mat,const MatType,Mat*); 56 extern PetscErrorCode MatDestroy_ML(Mat); 57 extern PetscErrorCode MatConvert_ML_SeqAIJ(ML_Operator*,Mat*); 58 extern PetscErrorCode MatConvert_ML_MPIAIJ(ML_Operator*,Mat*); 59 extern PetscErrorCode MatConvert_ML_SHELL(ML_Operator*,Mat*); 60 61 /* -------------------------------------------------------------------------- */ 62 /* 63 PCSetUp_ML - Prepares for the use of the ML preconditioner 64 by setting data structures and options. 65 66 Input Parameter: 67 . pc - the preconditioner context 68 69 Application Interface Routine: PCSetUp() 70 71 Notes: 72 The interface routine PCSetUp() is not usually called directly by 73 the user, but instead is called by PCApply() if necessary. 74 */ 75 76 #undef __FUNCT__ 77 #define __FUNCT__ "PCSetUp_ML" 78 static PetscErrorCode PCSetUp_ML(PC pc) 79 { 80 PetscErrorCode ierr; 81 PetscMPIInt size,rank; 82 FineGridCtx *PetscMLdata; 83 ML *ml_object; 84 ML_Aggregate *agg_object; 85 ML_Operator *mlmat; 86 PetscInt nlocal_allcols,Nlevels,mllevel,level,m,fine_level; 87 Mat A,Aloc; 88 GridCtx *gridctx; 89 PC_ML *pc_ml=PETSC_NULL; 90 PetscObjectContainer container; 91 92 PetscFunctionBegin; 93 ierr = PetscObjectQuery((PetscObject)pc,"PC_ML",(PetscObject *)&container);CHKERRQ(ierr); 94 if (container) { 95 ierr = PetscObjectContainerGetPointer(container,(void **)&pc_ml);CHKERRQ(ierr); 96 } else { 97 SETERRQ(PETSC_ERR_ARG_NULL,"Container does not exit"); 98 } 99 100 /* setup special features of PCML */ 101 /*--------------------------------*/ 102 /* covert A to Aloc to be used by ML at fine grid */ 103 A = pc->pmat; 104 ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr); 105 ierr = MPI_Comm_rank(A->comm,&rank);CHKERRQ(ierr); /* rm! */ 106 pc_ml->size = size; 107 if (size > 1){ 108 Aloc = PETSC_NULL; 109 ierr = MatConvert_MPIAIJ_ML(A,0,&Aloc);CHKERRQ(ierr); 110 } else { 111 Aloc = A; 112 } 113 114 /* create and initialize struct 'PetscMLdata' */ 115 ierr = PetscNew(FineGridCtx,&PetscMLdata);CHKERRQ(ierr); 116 PetscMLdata->A = A; 117 PetscMLdata->Aloc = Aloc; 118 ierr = PetscMalloc((Aloc->n+1)*sizeof(PetscScalar),&PetscMLdata->pwork);CHKERRQ(ierr); 119 pc_ml->PetscMLdata = PetscMLdata; 120 121 ierr = VecCreate(PETSC_COMM_SELF,&PetscMLdata->x);CHKERRQ(ierr); 122 if (size == 1){ 123 ierr = VecSetSizes(PetscMLdata->x,A->n,PETSC_DECIDE);CHKERRQ(ierr); 124 } else { 125 ierr = VecSetSizes(PetscMLdata->x,Aloc->n,PETSC_DECIDE);CHKERRQ(ierr); 126 } 127 ierr = VecSetType(PetscMLdata->x,VECSEQ);CHKERRQ(ierr); 128 129 ierr = VecCreate(PETSC_COMM_SELF,&PetscMLdata->y);CHKERRQ(ierr); 130 ierr = VecSetSizes(PetscMLdata->y,A->m,PETSC_DECIDE);CHKERRQ(ierr); 131 ierr = VecSetType(PetscMLdata->y,VECSEQ);CHKERRQ(ierr); 132 133 /* create ML discretization matrix at fine grid */ 134 ierr = MatGetSize(Aloc,&m,&nlocal_allcols);CHKERRQ(ierr); 135 ML_Create(&ml_object,pc_ml->MaxNlevels); 136 ML_Init_Amatrix(ml_object,0,m,m,PetscMLdata); 137 ML_Set_Amatrix_Getrow(ml_object,0,PetscML_getrow,PetscML_comm,nlocal_allcols); 138 ML_Set_Amatrix_Matvec(ml_object,0,PetscML_matvec); 139 140 /* aggregation */ 141 ML_Aggregate_Create(&agg_object); 142 ML_Aggregate_Set_MaxCoarseSize(agg_object,pc_ml->MaxCoarseSize); 143 /* set options */ 144 switch (pc_ml->CoarsenScheme) { 145 case 1: 146 ML_Aggregate_Set_CoarsenScheme_Coupled(agg_object);break; 147 case 2: 148 ML_Aggregate_Set_CoarsenScheme_MIS(agg_object);break; 149 case 3: 150 ML_Aggregate_Set_CoarsenScheme_METIS(agg_object);break; 151 } 152 ML_Aggregate_Set_Threshold(agg_object,pc_ml->Threshold); 153 ML_Aggregate_Set_DampingFactor(agg_object,pc_ml->DampingFactor); 154 if (pc_ml->SpectralNormScheme_Anorm){ 155 ML_Aggregate_Set_SpectralNormScheme_Anorm(agg_object); 156 } 157 158 Nlevels = ML_Gen_MGHierarchy_UsingAggregation(ml_object,0,ML_INCREASING,agg_object); 159 if (Nlevels<=0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Nlevels %d must > 0",Nlevels); 160 ierr = MGSetLevels(pc,Nlevels,PETSC_NULL);CHKERRQ(ierr); 161 pc_ml->ml_object = ml_object; 162 pc_ml->agg_object = agg_object; 163 164 ierr = PetscMalloc(Nlevels*sizeof(GridCtx),&gridctx);CHKERRQ(ierr); 165 fine_level = Nlevels - 1; 166 pc_ml->gridctx = gridctx; 167 pc_ml->fine_level = fine_level; 168 169 /* wrap ML matrices by PETSc shell matrices at coarsened grids. 170 Level 0 is the finest grid for ML, but coarsest for PETSc! */ 171 gridctx[fine_level].A = A; 172 level = fine_level - 1; 173 if (size == 1){ /* convert ML P, R and A into seqaij format */ 174 for (mllevel=1; mllevel<Nlevels; mllevel++){ 175 mlmat = &(ml_object->Pmat[mllevel]); 176 ierr = MatConvert_ML_SeqAIJ(mlmat,&gridctx[level].P);CHKERRQ(ierr); 177 mlmat = &(ml_object->Amat[mllevel]); 178 ierr = MatConvert_ML_SeqAIJ(mlmat,&gridctx[level].A);CHKERRQ(ierr); 179 mlmat = &(ml_object->Rmat[mllevel-1]); 180 ierr = MatConvert_ML_SeqAIJ(mlmat,&gridctx[level].R);CHKERRQ(ierr); 181 level--; 182 } 183 } else { /* convert ML P and R into shell format, ML A into mpiaij format */ 184 for (mllevel=1; mllevel<Nlevels; mllevel++){ 185 mlmat = &(ml_object->Pmat[mllevel]); 186 ierr = MatConvert_ML_SHELL(mlmat,&gridctx[level].P);CHKERRQ(ierr); 187 mlmat = &(ml_object->Rmat[mllevel-1]); 188 ierr = MatConvert_ML_SHELL(mlmat,&gridctx[level].R);CHKERRQ(ierr); 189 mlmat = &(ml_object->Amat[mllevel]); 190 ierr = MatConvert_ML_MPIAIJ(mlmat,&gridctx[level].A);CHKERRQ(ierr); 191 #ifdef TMP 192 /* 193 if (mllevel==1) { 194 ML_Operator_Print_UsingGlobalOrdering(mlmat,"Amat1",PETSC_NULL,PETSC_NULL); 195 } 196 */ 197 ierr = MatConvert_ML_SHELL(mlmat,&gridctx[level].A);CHKERRQ(ierr); 198 if (mllevel>0){ 199 PetscMLdata->mlmat = &(ml_object->Amat[mllevel]); 200 Mat A_tmp; 201 ierr = MatConvert_ML_MPIAIJ(PetscMLdata,&A_tmp);CHKERRQ(ierr); 202 /* ierr = MatView(A_tmp,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); */ 203 204 Vec x,yy1,yy2; 205 PetscInt am,an; 206 ierr = MatGetLocalSize(A_tmp,&am,&an);CHKERRQ(ierr); 207 ierr = VecCreate(PETSC_COMM_WORLD,&x);CHKERRQ(ierr); 208 ierr = VecSetSizes(x,an,PETSC_DECIDE);CHKERRQ(ierr); 209 ierr = VecSetFromOptions(x);CHKERRQ(ierr); 210 ierr = VecCreate(PETSC_COMM_WORLD,&yy1);CHKERRQ(ierr); 211 ierr = VecSetSizes(yy1,am,PETSC_DECIDE);CHKERRQ(ierr); 212 ierr = VecSetFromOptions(yy1);CHKERRQ(ierr); 213 ierr = VecDuplicate(yy1,&yy2);CHKERRQ(ierr); 214 215 PetscRandom rd; 216 PetscReal rnorm; 217 PetscScalar mone = -1.0; 218 ierr = PetscRandomCreate(PETSC_COMM_WORLD,RANDOM_DEFAULT,&rd);CHKERRQ(ierr); 219 ierr = VecSetRandom(rd,x);CHKERRQ(ierr); 220 ierr = MatMult(gridctx[level].A,x,yy1);CHKERRQ(ierr); 221 ierr = MatMult(A_tmp,x,yy2);CHKERRQ(ierr); 222 ierr = VecAXPY(&mone,yy1,yy2);CHKERRQ(ierr); 223 ierr = VecNorm(yy2,NORM_2,&rnorm);CHKERRQ(ierr); 224 printf(" [%d] mllevel: %d rnorm %g\n",rank,mllevel,rnorm); 225 226 ierr = MatDestroy(A_tmp);CHKERRQ(ierr); 227 ierr = VecDestroy(x);CHKERRQ(ierr); 228 ierr = VecDestroy(yy1);CHKERRQ(ierr); 229 ierr = VecDestroy(yy2);CHKERRQ(ierr); 230 ierr = PetscRandomDestroy(rd);CHKERRQ(ierr); 231 } 232 #endif 233 level--; 234 } 235 } 236 237 /* create coarse level and the interpolation between the levels */ 238 level = fine_level; 239 while ( level >= 0 ){ 240 if (level != fine_level){ 241 ierr = VecCreate(gridctx[level].A->comm,&gridctx[level].x);CHKERRQ(ierr); 242 ierr = VecSetSizes(gridctx[level].x,gridctx[level].A->n,PETSC_DECIDE);CHKERRQ(ierr); 243 ierr = VecSetType(gridctx[level].x,VECMPI);CHKERRQ(ierr); 244 ierr = MGSetX(pc,level,gridctx[level].x);CHKERRQ(ierr); 245 246 ierr = VecCreate(gridctx[level].A->comm,&gridctx[level].b);CHKERRQ(ierr); 247 ierr = VecSetSizes(gridctx[level].b,gridctx[level].A->m,PETSC_DECIDE);CHKERRQ(ierr); 248 ierr = VecSetType(gridctx[level].b,VECMPI);CHKERRQ(ierr); 249 ierr = MGSetRhs(pc,level,gridctx[level].b);CHKERRQ(ierr); 250 } 251 ierr = VecCreate(gridctx[level].A->comm,&gridctx[level].r);CHKERRQ(ierr); 252 ierr = VecSetSizes(gridctx[level].r,gridctx[level].A->m,PETSC_DECIDE);CHKERRQ(ierr); 253 ierr = VecSetType(gridctx[level].r,VECMPI);CHKERRQ(ierr); 254 ierr = MGSetR(pc,level,gridctx[level].r);CHKERRQ(ierr); 255 256 if (level == 0){ 257 ierr = MGGetCoarseSolve(pc,&gridctx[level].ksp);CHKERRQ(ierr); 258 } else { 259 ierr = MGGetSmoother(pc,level,&gridctx[level].ksp);CHKERRQ(ierr); 260 ierr = MGSetResidual(pc,level,MGDefaultResidual,gridctx[level].A);CHKERRQ(ierr); 261 if (level == fine_level){ 262 ierr = KSPSetOptionsPrefix(gridctx[level].ksp,"mg_fine_");CHKERRQ(ierr); 263 ierr = MGSetR(pc,level,gridctx[level].r);CHKERRQ(ierr); 264 } 265 } 266 ierr = KSPSetOperators(gridctx[level].ksp,gridctx[level].A,gridctx[level].A,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 267 268 if (level < fine_level){ 269 ierr = MGSetInterpolate(pc,level+1,gridctx[level].P);CHKERRQ(ierr); 270 ierr = MGSetRestriction(pc,level+1,gridctx[level].R);CHKERRQ(ierr); 271 } 272 level--; 273 } 274 275 /* now call PCSetUp_MG() */ 276 /*--------------------------------*/ 277 ierr = (*pc_ml->PCSetUp)(pc);CHKERRQ(ierr); 278 PetscFunctionReturn(0); 279 } 280 281 #undef __FUNCT__ 282 #define __FUNCT__ "PetscObjectContainerDestroy_PC_ML" 283 PetscErrorCode PetscObjectContainerDestroy_PC_ML(void *ptr) 284 { 285 PetscErrorCode ierr; 286 PC_ML *pc_ml = (PC_ML*)ptr; 287 PetscInt level; 288 289 PetscFunctionBegin; 290 if (pc_ml->size > 1){ierr = MatDestroy(pc_ml->PetscMLdata->Aloc);CHKERRQ(ierr);} 291 ML_Aggregate_Destroy(&pc_ml->agg_object); 292 ML_Destroy(&pc_ml->ml_object); 293 294 ierr = PetscFree(pc_ml->PetscMLdata->pwork);CHKERRQ(ierr); 295 ierr = VecDestroy(pc_ml->PetscMLdata->x);CHKERRQ(ierr); 296 ierr = VecDestroy(pc_ml->PetscMLdata->y);CHKERRQ(ierr); 297 ierr = PetscFree(pc_ml->PetscMLdata);CHKERRQ(ierr); 298 299 level = pc_ml->fine_level; 300 while ( level>= 0 ){ 301 if (level != pc_ml->fine_level){ 302 ierr = MatDestroy(pc_ml->gridctx[level].A);CHKERRQ(ierr); 303 ierr = MatDestroy(pc_ml->gridctx[level].P);CHKERRQ(ierr); 304 ierr = MatDestroy(pc_ml->gridctx[level].R);CHKERRQ(ierr); 305 ierr = VecDestroy(pc_ml->gridctx[level].x);CHKERRQ(ierr); 306 ierr = VecDestroy(pc_ml->gridctx[level].b);CHKERRQ(ierr); 307 } 308 ierr = VecDestroy(pc_ml->gridctx[level].r);CHKERRQ(ierr); 309 level--; 310 } 311 ierr = PetscFree(pc_ml->gridctx);CHKERRQ(ierr); 312 ierr = PetscFree(pc_ml);CHKERRQ(ierr); 313 PetscFunctionReturn(0); 314 } 315 /* -------------------------------------------------------------------------- */ 316 /* 317 PCDestroy_ML - Destroys the private context for the ML preconditioner 318 that was created with PCCreate_ML(). 319 320 Input Parameter: 321 . pc - the preconditioner context 322 323 Application Interface Routine: PCDestroy() 324 */ 325 #undef __FUNCT__ 326 #define __FUNCT__ "PCDestroy_ML" 327 static PetscErrorCode PCDestroy_ML(PC pc) 328 { 329 PetscErrorCode ierr; 330 PC_ML *pc_ml=PETSC_NULL; 331 PetscObjectContainer container; 332 333 PetscFunctionBegin; 334 ierr = PetscObjectQuery((PetscObject)pc,"PC_ML",(PetscObject *)&container);CHKERRQ(ierr); 335 if (container) { 336 ierr = PetscObjectContainerGetPointer(container,(void **)&pc_ml);CHKERRQ(ierr); 337 pc->ops->destroy = pc_ml->PCDestroy; 338 } else { 339 SETERRQ(PETSC_ERR_ARG_NULL,"Container does not exit"); 340 } 341 /* detach pc and PC_ML and dereference container */ 342 ierr = PetscObjectCompose((PetscObject)pc,"PC_ML",0);CHKERRQ(ierr); 343 ierr = (*pc->ops->destroy)(pc);CHKERRQ(ierr); 344 345 ierr = PetscObjectContainerDestroy(container);CHKERRQ(ierr); 346 PetscFunctionReturn(0); 347 } 348 349 #undef __FUNCT__ 350 #define __FUNCT__ "PCSetFromOptions_ML" 351 static PetscErrorCode PCSetFromOptions_ML(PC pc) 352 { 353 PetscErrorCode ierr; 354 PetscInt indx,m,PrintLevel,MaxNlevels,MaxCoarseSize; 355 PetscReal Threshold,DampingFactor; 356 PetscTruth flg; 357 const char *type[] = {"additive","multiplicative","full","cascade","kascade"}; 358 const char *scheme[] = {"Uncoupled","Coupled","MIS","METIS"}; 359 PC_ML *pc_ml=PETSC_NULL; 360 PetscObjectContainer container; 361 362 PetscFunctionBegin; 363 ierr = PetscObjectQuery((PetscObject)pc,"PC_ML",(PetscObject *)&container);CHKERRQ(ierr); 364 if (container) { 365 ierr = PetscObjectContainerGetPointer(container,(void **)&pc_ml);CHKERRQ(ierr); 366 } else { 367 SETERRQ(PETSC_ERR_ARG_NULL,"Container does not exit"); 368 } 369 ierr = PetscOptionsHead("MG options");CHKERRQ(ierr); 370 /* inherited MG options */ 371 ierr = PetscOptionsInt("-pc_mg_cycles","1 for V cycle, 2 for W-cycle","MGSetCycles",1,&m,&flg);CHKERRQ(ierr); 372 if (flg) { 373 ierr = MGSetCycles(pc,m);CHKERRQ(ierr); 374 } 375 ierr = PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","MGSetNumberSmoothUp",1,&m,&flg);CHKERRQ(ierr); 376 if (flg) { 377 ierr = MGSetNumberSmoothUp(pc,m);CHKERRQ(ierr); 378 } 379 ierr = PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","MGSetNumberSmoothDown",1,&m,&flg);CHKERRQ(ierr); 380 if (flg) { 381 ierr = MGSetNumberSmoothDown(pc,m);CHKERRQ(ierr); 382 } 383 ierr = PetscOptionsEList("-pc_mg_type","Multigrid type","MGSetType",type,5,type[1],&indx,&flg);CHKERRQ(ierr); 384 if (flg) { 385 MGType mg = (MGType) 0; 386 switch (indx) { 387 case 0: 388 mg = MGADDITIVE; 389 break; 390 case 1: 391 mg = MGMULTIPLICATIVE; 392 break; 393 case 2: 394 mg = MGFULL; 395 break; 396 case 3: 397 mg = MGKASKADE; 398 break; 399 case 4: 400 mg = MGKASKADE; 401 break; 402 } 403 ierr = MGSetType(pc,mg);CHKERRQ(ierr); 404 } 405 ierr = PetscOptionsTail();CHKERRQ(ierr); 406 407 /* ML options */ 408 ierr = PetscOptionsHead("ML options");CHKERRQ(ierr); 409 /* set defaults */ 410 PrintLevel = 0; 411 MaxNlevels = 10; 412 MaxCoarseSize = 1; 413 indx = 0; 414 Threshold = 0.0; 415 DampingFactor = 4.0/3.0; 416 417 ierr = PetscOptionsInt("-pc_ml_PrintLevel","Print level","ML_Set_PrintLevel",PrintLevel,&PrintLevel,PETSC_NULL);CHKERRQ(ierr); 418 ML_Set_PrintLevel(PrintLevel); 419 420 ierr = PetscOptionsInt("-pc_ml_maxNlevels","Maximum number of levels","None",MaxNlevels,&MaxNlevels,PETSC_NULL);CHKERRQ(ierr); 421 pc_ml->MaxNlevels = MaxNlevels; 422 423 ierr = PetscOptionsInt("-pc_ml_maxCoarseSize","Maximum coarsest mesh size","ML_Aggregate_Set_MaxCoarseSize",MaxCoarseSize,&MaxCoarseSize,PETSC_NULL);CHKERRQ(ierr); 424 pc_ml->MaxCoarseSize = MaxCoarseSize; 425 426 ierr = PetscOptionsEList("-pc_ml_CoarsenScheme","Aggregate Coarsen Scheme","ML_Aggregate_Set_CoarsenScheme_*",scheme,4,scheme[0],&indx,PETSC_NULL);CHKERRQ(ierr); 427 pc_ml->CoarsenScheme = indx; 428 429 ierr = PetscOptionsReal("-pc_ml_DampingFactor","P damping factor","ML_Aggregate_Set_DampingFactor",DampingFactor,&DampingFactor,PETSC_NULL);CHKERRQ(ierr); 430 pc_ml->DampingFactor = DampingFactor; 431 432 ierr = PetscOptionsReal("-pc_ml_Threshold","Smoother drop tol","ML_Aggregate_Set_Threshold",Threshold,&Threshold,PETSC_NULL);CHKERRQ(ierr); 433 pc_ml->Threshold = Threshold; 434 435 ierr = PetscOptionsLogical("-pc_ml_SpectralNormScheme_Anorm","Method used for estimating spectral radius","ML_Aggregate_Set_SpectralNormScheme_Anorm",PETSC_FALSE,&pc_ml->SpectralNormScheme_Anorm,PETSC_FALSE); 436 437 ierr = PetscOptionsTail();CHKERRQ(ierr); 438 PetscFunctionReturn(0); 439 } 440 441 /* -------------------------------------------------------------------------- */ 442 /* 443 PCCreate_ML - Creates a ML preconditioner context, PC_ML, 444 and sets this as the private data within the generic preconditioning 445 context, PC, that was created within PCCreate(). 446 447 Input Parameter: 448 . pc - the preconditioner context 449 450 Application Interface Routine: PCCreate() 451 */ 452 453 /*MC 454 PCML - Use geometric multigrid preconditioning. This preconditioner requires you provide 455 fine grid discretization matrix. The coarser grid matrices and restriction/interpolation 456 operators are computed by ML and wrapped as PETSc shell matrices. 457 458 Options Database Key: (not done yet!) 459 + -pc_mg_maxlevels <nlevels> - maximum number of levels including finest 460 . -pc_mg_cycles 1 or 2 - for V or W-cycle 461 . -pc_mg_smoothup <n> - number of smoothing steps after interpolation 462 . -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator 463 . -pc_mg_type <additive,multiplicative,full,cascade> - multiplicative is the default 464 . -pc_mg_monitor - print information on the multigrid convergence 465 - -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices 466 to the Socket viewer for reading from Matlab. 467 468 Level: intermediate 469 470 Concepts: multigrid 471 472 .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, 473 MGSetLevels(), MGGetLevels(), MGSetType(), MPSetCycles(), MGSetNumberSmoothDown(), 474 MGSetNumberSmoothUp(), MGGetCoarseSolve(), MGSetResidual(), MGSetInterpolation(), 475 MGSetRestriction(), MGGetSmoother(), MGGetSmootherUp(), MGGetSmootherDown(), 476 MGSetCyclesOnLevel(), MGSetRhs(), MGSetX(), MGSetR() 477 M*/ 478 479 EXTERN_C_BEGIN 480 #undef __FUNCT__ 481 #define __FUNCT__ "PCCreate_ML" 482 PetscErrorCode PCCreate_ML(PC pc) 483 { 484 PetscErrorCode ierr; 485 PC_ML *pc_ml; 486 PetscObjectContainer container; 487 488 PetscFunctionBegin; 489 /* initialize pc as PCMG */ 490 ierr = PCSetType(pc,PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */ 491 492 /* create a supporting struct and attach it to pc */ 493 ierr = PetscNew(PC_ML,&pc_ml);CHKERRQ(ierr); 494 ierr = PetscObjectContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr); 495 ierr = PetscObjectContainerSetPointer(container,pc_ml);CHKERRQ(ierr); 496 ierr = PetscObjectContainerSetUserDestroy(container,PetscObjectContainerDestroy_PC_ML);CHKERRQ(ierr); 497 ierr = PetscObjectCompose((PetscObject)pc,"PC_ML",(PetscObject)container);CHKERRQ(ierr); 498 499 pc_ml->PCSetUp = pc->ops->setup; 500 pc_ml->PCDestroy = pc->ops->destroy; 501 502 /* overwrite the pointers of PCMG by the functions of PCML */ 503 pc->ops->setfromoptions = PCSetFromOptions_ML; 504 pc->ops->setup = PCSetUp_ML; 505 pc->ops->destroy = PCDestroy_ML; 506 PetscFunctionReturn(0); 507 } 508 EXTERN_C_END 509 510 int PetscML_getrow(void *ML_data, int N_requested_rows, int requested_rows[], 511 int allocated_space, int columns[], double values[], int row_lengths[]) 512 { 513 PetscErrorCode ierr; 514 Mat Aloc; 515 Mat_SeqAIJ *a; 516 PetscInt m,i,j,k=0,row,*aj; 517 PetscScalar *aa; 518 FineGridCtx *ml=(FineGridCtx*)ML_data; 519 520 Aloc = ml->Aloc; 521 a = (Mat_SeqAIJ*)Aloc->data; 522 ierr = MatGetSize(Aloc,&m,PETSC_NULL);CHKERRQ(ierr); 523 524 for (i = 0; i<N_requested_rows; i++) { 525 row = requested_rows[i]; 526 row_lengths[i] = a->ilen[row]; 527 if (allocated_space < k+row_lengths[i]) return(0); 528 if ( (row >= 0) || (row <= (m-1)) ) { 529 aj = a->j + a->i[row]; 530 aa = a->a + a->i[row]; 531 for (j=0; j<row_lengths[i]; j++){ 532 columns[k] = aj[j]; 533 values[k++] = aa[j]; 534 } 535 } 536 } 537 return(1); 538 } 539 540 int PetscML_matvec(void *ML_data,int in_length,double p[],int out_length,double ap[]) 541 { 542 PetscErrorCode ierr; 543 FineGridCtx *ml=(FineGridCtx*)ML_data; 544 Mat A=ml->A, Aloc=ml->Aloc; 545 PetscMPIInt size; 546 PetscScalar *pwork=ml->pwork; 547 PetscInt i; 548 549 ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr); 550 if (size == 1){ 551 ierr = VecPlaceArray(ml->x,p);CHKERRQ(ierr); 552 } else { 553 for (i=0; i<in_length; i++) pwork[i] = p[i]; 554 PetscML_comm(pwork,ml); 555 ierr = VecPlaceArray(ml->x,pwork);CHKERRQ(ierr); 556 } 557 ierr = VecPlaceArray(ml->y,ap);CHKERRQ(ierr); 558 ierr = MatMult(Aloc,ml->x,ml->y);CHKERRQ(ierr); 559 return 0; 560 } 561 562 int PetscML_comm(double p[],void *ML_data) 563 { 564 PetscErrorCode ierr; 565 FineGridCtx *ml=(FineGridCtx*)ML_data; 566 Mat A=ml->A; 567 Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data; 568 PetscMPIInt size; 569 PetscInt i,in_length=A->m,out_length=ml->Aloc->n; 570 PetscScalar *array; 571 572 ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr); 573 if (size == 1) return 0; 574 575 ierr = VecPlaceArray(ml->y,p);CHKERRQ(ierr); 576 ierr = VecScatterBegin(ml->y,a->lvec,INSERT_VALUES,SCATTER_FORWARD,a->Mvctx);CHKERRQ(ierr); 577 ierr = VecScatterEnd(ml->y,a->lvec,INSERT_VALUES,SCATTER_FORWARD,a->Mvctx);CHKERRQ(ierr); 578 ierr = VecGetArray(a->lvec,&array);CHKERRQ(ierr); 579 for (i=in_length; i<out_length; i++){ 580 p[i] = array[i-in_length]; 581 } 582 return 0; 583 } 584 #undef __FUNCT__ 585 #define __FUNCT__ "MatMult_ML" 586 PetscErrorCode MatMult_ML(Mat A,Vec x,Vec y) 587 { 588 PetscErrorCode ierr; 589 Mat_MLShell *shell; 590 PetscScalar *xarray,*yarray; 591 PetscInt x_length,y_length; 592 593 PetscFunctionBegin; 594 ierr = MatShellGetContext(A,(void *)&shell);CHKERRQ(ierr); 595 ierr = VecGetArray(x,&xarray);CHKERRQ(ierr); 596 ierr = VecGetArray(y,&yarray);CHKERRQ(ierr); 597 x_length = shell->mlmat->invec_leng; 598 y_length = shell->mlmat->outvec_leng; 599 600 ML_Operator_Apply(shell->mlmat,x_length,xarray,y_length,yarray); 601 602 ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr); 603 ierr = VecRestoreArray(y,&yarray);CHKERRQ(ierr); 604 PetscFunctionReturn(0); 605 } 606 /* MatMultAdd_ML - Compute y = w + A*x */ 607 #undef __FUNCT__ 608 #define __FUNCT__ "MatMultAdd_ML" 609 PetscErrorCode MatMultAdd_ML(Mat A,Vec x,Vec w,Vec y) 610 { 611 PetscErrorCode ierr; 612 Mat_MLShell *shell; 613 PetscScalar *xarray,*yarray; 614 const PetscScalar one=1.0; 615 PetscInt x_length,y_length; 616 617 PetscFunctionBegin; 618 ierr = MatShellGetContext(A,(void *)&shell);CHKERRQ(ierr); 619 ierr = VecGetArray(x,&xarray);CHKERRQ(ierr); 620 ierr = VecGetArray(y,&yarray);CHKERRQ(ierr); 621 622 x_length = shell->mlmat->invec_leng; 623 y_length = shell->mlmat->outvec_leng; 624 625 ML_Operator_Apply(shell->mlmat,x_length,xarray,y_length,yarray); 626 627 ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr); 628 ierr = VecRestoreArray(y,&yarray);CHKERRQ(ierr); 629 ierr = VecAXPY(&one,w,y);CHKERRQ(ierr); 630 631 PetscFunctionReturn(0); 632 } 633 634 /* newtype is ignored because "ml" is not listed under Petsc MatType yet */ 635 #undef __FUNCT__ 636 #define __FUNCT__ "MatConvert_MPIAIJ_ML" 637 PetscErrorCode MatConvert_MPIAIJ_ML(Mat A,const MatType newtype,Mat *Aloc) 638 { 639 PetscErrorCode ierr; 640 Mat_MPIAIJ *mpimat=(Mat_MPIAIJ*)A->data; 641 Mat_SeqAIJ *mat,*a=(Mat_SeqAIJ*)(mpimat->A)->data,*b=(Mat_SeqAIJ*)(mpimat->B)->data; 642 PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j; 643 PetscScalar *aa=a->a,*ba=b->a,*ca; 644 PetscInt am=A->m,an=A->n,i,j,k; 645 PetscInt *ci,*cj,ncols; 646 MatReuse scall=MAT_INITIAL_MATRIX; 647 648 PetscFunctionBegin; 649 if (am != an) SETERRQ2(PETSC_ERR_ARG_WRONG,"A must have a square diagonal portion, am: %d != an: %d",am,an); 650 651 if (*Aloc) scall = MAT_REUSE_MATRIX; 652 if (scall == MAT_INITIAL_MATRIX){ 653 ierr = PetscMalloc((1+am)*sizeof(PetscInt),&ci);CHKERRQ(ierr); 654 ci[0] = 0; 655 for (i=0; i<am; i++){ 656 ci[i+1] = ci[i] + (ai[i+1] - ai[i]) + (bi[i+1] - bi[i]); 657 } 658 ierr = PetscMalloc((1+ci[am])*sizeof(PetscInt),&cj);CHKERRQ(ierr); 659 ierr = PetscMalloc((1+ci[am])*sizeof(PetscScalar),&ca);CHKERRQ(ierr); 660 661 k = 0; 662 for (i=0; i<am; i++){ 663 /* diagonal portion of A */ 664 ncols = ai[i+1] - ai[i]; 665 for (j=0; j<ncols; j++) { 666 cj[k] = *aj++; 667 ca[k++] = *aa++; 668 } 669 /* off-diagonal portion of A */ 670 ncols = bi[i+1] - bi[i]; 671 for (j=0; j<ncols; j++) { 672 cj[k] = an + (*bj); bj++; 673 ca[k++] = *ba++; 674 } 675 } 676 if (k != ci[am]) SETERRQ2(PETSC_ERR_ARG_WRONG,"k: %d != ci[am]: %d",k,ci[am]); 677 678 /* put together the new matrix */ 679 an = mpimat->A->n+mpimat->B->n; 680 ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,an,ci,cj,ca,Aloc);CHKERRQ(ierr); 681 682 /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 683 /* Since these are PETSc arrays, change flags to free them as necessary. */ 684 mat = (Mat_SeqAIJ*)(*Aloc)->data; 685 mat->freedata = PETSC_TRUE; 686 mat->nonew = 0; 687 } else if (scall == MAT_REUSE_MATRIX){ 688 mat=(Mat_SeqAIJ*)(*Aloc)->data; 689 ci = mat->i; cj = mat->j; ca = mat->a; 690 for (i=0; i<am; i++) { 691 /* diagonal portion of A */ 692 ncols = ai[i+1] - ai[i]; 693 for (j=0; j<ncols; j++) *ca++ = *aa++; 694 /* off-diagonal portion of A */ 695 ncols = bi[i+1] - bi[i]; 696 for (j=0; j<ncols; j++) *ca++ = *ba++; 697 } 698 } else { 699 SETERRQ1(PETSC_ERR_ARG_WRONG,"Invalid MatReuse %d",(int)scall); 700 } 701 PetscFunctionReturn(0); 702 } 703 extern PetscErrorCode MatDestroy_Shell(Mat); 704 #undef __FUNCT__ 705 #define __FUNCT__ "MatDestroy_ML" 706 PetscErrorCode MatDestroy_ML(Mat A) 707 { 708 PetscErrorCode ierr; 709 Mat_MLShell *shell; 710 711 PetscFunctionBegin; 712 ierr = MatShellGetContext(A,(void *)&shell);CHKERRQ(ierr); 713 ierr = VecDestroy(shell->y);CHKERRQ(ierr); 714 ierr = PetscFree(shell);CHKERRQ(ierr); 715 ierr = MatDestroy_Shell(A);CHKERRQ(ierr); 716 PetscFunctionReturn(0); 717 } 718 719 extern PetscErrorCode PetscSortIntWithScalarArray(PetscInt,PetscInt [],PetscScalar []); 720 #undef __FUNCT__ 721 #define __FUNCT__ "MatConvert_ML_SeqAIJ" 722 PetscErrorCode MatConvert_ML_SeqAIJ(ML_Operator *mlmat,Mat *newmat) 723 { 724 struct ML_CSR_MSRdata *matdata = (struct ML_CSR_MSRdata *)mlmat->data; 725 PetscErrorCode ierr; 726 PetscInt m=mlmat->outvec_leng,n=mlmat->invec_leng,nnz[m],nz_max; 727 PetscInt *ml_cols=matdata->columns,*aj,i,j,k; 728 PetscScalar *ml_vals=matdata->values,*aa; 729 730 PetscFunctionBegin; 731 if ( mlmat->getrow == NULL) SETERRQ(PETSC_ERR_ARG_NULL,"mlmat->getrow = NULL"); 732 if (m != n){ /* ML Pmat and Rmat are in CSR format. Pass array pointers into SeqAIJ matrix */ 733 ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,m,n,matdata->rowptr,ml_cols,ml_vals,newmat);CHKERRQ(ierr); 734 PetscFunctionReturn(0); 735 } 736 737 /* ML Amat is in MSR format. Copy its data into SeqAIJ matrix */ 738 ierr = MatCreate(PETSC_COMM_SELF,m,n,PETSC_DECIDE,PETSC_DECIDE,newmat);CHKERRQ(ierr); 739 ierr = MatSetType(*newmat,MATSEQAIJ);CHKERRQ(ierr); 740 741 nz_max = 0; 742 for (i=0; i<m; i++) { 743 nnz[i] = ml_cols[i+1] - ml_cols[i] + 1; 744 if (nnz[i] > nz_max) nz_max = nnz[i]; 745 } 746 ierr = MatSeqAIJSetPreallocation(*newmat,0,nnz);CHKERRQ(ierr); 747 ierr = MatSetOption(*newmat,MAT_COLUMNS_SORTED);CHKERRQ(ierr); /* check! */ 748 749 nz_max++; 750 ierr = PetscMalloc(nz_max*(sizeof(int)+sizeof(PetscScalar)),&aj);CHKERRQ(ierr); 751 aa = (PetscScalar*)(aj + nz_max); 752 753 for (i=0; i<m; i++){ 754 k = 0; 755 /* diagonal entry */ 756 aj[k] = i; aa[k++] = ml_vals[i]; 757 /* off diagonal entries */ 758 for (j=ml_cols[i]; j<ml_cols[i+1]; j++){ 759 aj[k] = ml_cols[j]; aa[k++] = ml_vals[j]; 760 } 761 /* sort aj and aa */ 762 ierr = PetscSortIntWithScalarArray(nnz[i],aj,aa);CHKERRQ(ierr); 763 ierr = MatSetValues(*newmat,1,&i,nnz[i],aj,aa,INSERT_VALUES);CHKERRQ(ierr); 764 } 765 ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 766 ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 767 ierr = PetscFree(aj);CHKERRQ(ierr); 768 PetscFunctionReturn(0); 769 } 770 771 #undef __FUNCT__ 772 #define __FUNCT__ "MatConvert_ML_SHELL" 773 PetscErrorCode MatConvert_ML_SHELL(ML_Operator *mlmat,Mat *newmat) 774 { 775 PetscErrorCode ierr; 776 PetscInt m,n; 777 ML_Comm *MLcomm; 778 Mat_MLShell *shellctx; 779 780 PetscFunctionBegin; 781 m = mlmat->outvec_leng; 782 n = mlmat->invec_leng; 783 if (!m || !n){ 784 newmat = PETSC_NULL; 785 } else { 786 MLcomm = mlmat->comm; 787 ierr = PetscNew(Mat_MLShell,&shellctx);CHKERRQ(ierr); 788 ierr = MatCreateShell(MLcomm->USR_comm,m,n,PETSC_DETERMINE,PETSC_DETERMINE,shellctx,newmat);CHKERRQ(ierr); 789 ierr = MatShellSetOperation(*newmat,MATOP_MULT,(void *)MatMult_ML);CHKERRQ(ierr); 790 ierr = MatShellSetOperation(*newmat,MATOP_MULT_ADD,(void *)MatMultAdd_ML);CHKERRQ(ierr); 791 shellctx->A = *newmat; 792 shellctx->mlmat = mlmat; 793 ierr = VecCreate(PETSC_COMM_WORLD,&shellctx->y);CHKERRQ(ierr); 794 ierr = VecSetSizes(shellctx->y,m,PETSC_DECIDE);CHKERRQ(ierr); 795 ierr = VecSetFromOptions(shellctx->y);CHKERRQ(ierr); 796 (*newmat)->ops->destroy = MatDestroy_ML; 797 } 798 PetscFunctionReturn(0); 799 } 800 801 #undef __FUNCT__ 802 #define __FUNCT__ "MatConvert_ML_MPIAIJ" 803 PetscErrorCode MatConvert_ML_MPIAIJ(ML_Operator *mlmat,Mat *newmat) 804 { 805 struct ML_CSR_MSRdata *matdata = (struct ML_CSR_MSRdata *)mlmat->data; 806 PetscInt *ml_cols=matdata->columns,*aj; 807 PetscScalar *ml_vals=matdata->values,*aa; 808 PetscErrorCode ierr; 809 PetscInt i,j,k,*gordering; 810 PetscInt m=mlmat->outvec_leng,n,nnzA[m],nnzB[m],nnz[m],nz_max,row; 811 Mat A; 812 813 PetscFunctionBegin; 814 if (mlmat->getrow == NULL) SETERRQ(PETSC_ERR_ARG_NULL,"mlmat->getrow = NULL"); 815 n = mlmat->invec_leng; 816 if (m != n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"m %d must equal to n %d",m,n); 817 818 ierr = MatCreate(mlmat->comm->USR_comm,m,n,PETSC_DECIDE,PETSC_DECIDE,&A);CHKERRQ(ierr); 819 ierr = MatSetType(A,MATMPIAIJ);CHKERRQ(ierr); 820 nz_max = 0; 821 for (i=0; i<m; i++){ 822 nnz[i] = ml_cols[i+1] - ml_cols[i] + 1; 823 if (nz_max < nnz[i]) nz_max = nnz[i]; 824 nnzA[i] = 1; /* diag */ 825 for (j=ml_cols[i]; j<ml_cols[i+1]; j++){ 826 if (ml_cols[j] < m) nnzA[i]++; 827 } 828 nnzB[i] = nnz[i] - nnzA[i]; 829 } 830 ierr = MatMPIAIJSetPreallocation(A,0,nnzA,0,nnzB);CHKERRQ(ierr); 831 832 /* insert mat values -- remap row and column indices */ 833 nz_max++; 834 ierr = PetscMalloc(nz_max*(sizeof(int)+sizeof(PetscScalar)),&aj);CHKERRQ(ierr); 835 aa = (PetscScalar*)(aj + nz_max); 836 ML_build_global_numbering(mlmat,mlmat->comm,&gordering); 837 for (i=0; i<m; i++){ 838 row = gordering[i]; 839 k = 0; 840 /* diagonal entry */ 841 aj[k] = row; aa[k++] = ml_vals[i]; 842 /* off diagonal entries */ 843 for (j=ml_cols[i]; j<ml_cols[i+1]; j++){ 844 aj[k] = gordering[ml_cols[j]]; aa[k++] = ml_vals[j]; 845 } 846 ierr = MatSetValues(A,1,&row,nnz[i],aj,aa,INSERT_VALUES);CHKERRQ(ierr); 847 } 848 ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 849 ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 850 *newmat = A; 851 852 ierr = PetscFree(aj);CHKERRQ(ierr); 853 PetscFunctionReturn(0); 854 } 855