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