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