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