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