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