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