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