15582bec1SHong Zhang /* 25582bec1SHong Zhang Provides an interface to the ML 3.0 smoothed Aggregation 35582bec1SHong Zhang */ 45582bec1SHong Zhang #include "src/ksp/pc/pcimpl.h" /*I "petscpc.h" I*/ 55582bec1SHong Zhang #include "src/ksp/pc/impls/mg/mgimpl.h" /*I "petscmg.h" I*/ 65582bec1SHong Zhang #include "src/mat/impls/aij/seq/aij.h" 75582bec1SHong Zhang #include "src/mat/impls/aij/mpi/mpiaij.h" 85582bec1SHong Zhang EXTERN_C_BEGIN 95582bec1SHong Zhang #include <math.h> 105582bec1SHong Zhang #include "ml_include.h" 115582bec1SHong Zhang EXTERN_C_END 125582bec1SHong Zhang 135582bec1SHong Zhang /* The context (data structure) at each grid level */ 145582bec1SHong Zhang typedef struct { 155582bec1SHong Zhang Vec x,b,r; /* global vectors */ 165582bec1SHong Zhang Mat A,P,R; 175582bec1SHong Zhang KSP ksp; 185582bec1SHong Zhang } GridCtx; 195582bec1SHong Zhang 205582bec1SHong Zhang /* The context used to input PETSc matrix into ML at fine grid */ 215582bec1SHong Zhang typedef struct { 225582bec1SHong Zhang Mat A,Aloc; 235582bec1SHong Zhang ML_Operator *mlmat; 245582bec1SHong Zhang PetscScalar *pwork; /* tmp array used by PetscML_comm() */ 255582bec1SHong Zhang PetscInt rlen_max,*cols; /* used by MatConvert_ML_SeqAIJ() */ 265582bec1SHong Zhang PetscScalar *vals; /* used by MatConvert_ML_SeqAIJ() */ 275582bec1SHong Zhang } FineGridCtx; 285582bec1SHong Zhang 295582bec1SHong Zhang /* The context associates a ML matrix with a PETSc shell matrix */ 305582bec1SHong Zhang typedef struct { 315582bec1SHong Zhang Mat A; /* PETSc shell matrix associated with mlmat */ 325582bec1SHong Zhang ML_Operator *mlmat; /* ML matrix assorciated with A */ 335582bec1SHong Zhang Vec y; 345582bec1SHong Zhang } Mat_MLShell; 355582bec1SHong Zhang 365582bec1SHong Zhang /* Private context for the ML preconditioner */ 375582bec1SHong Zhang typedef struct { 385582bec1SHong Zhang ML *ml_object; 395582bec1SHong Zhang ML_Aggregate *agg_object; 405582bec1SHong Zhang GridCtx *gridctx; 415582bec1SHong Zhang FineGridCtx *PetscMLdata; 425582bec1SHong Zhang PetscInt fine_level,MaxNlevels,MaxCoarseSize,CoarsenScheme; 435582bec1SHong Zhang PetscReal Threshold,DampingFactor; 445582bec1SHong Zhang PetscTruth SpectralNormScheme_Anorm; 455582bec1SHong Zhang PetscMPIInt size; 465582bec1SHong Zhang 475582bec1SHong Zhang PetscErrorCode (*PCSetUp)(PC); 485582bec1SHong Zhang PetscErrorCode (*PCDestroy)(PC); 495582bec1SHong Zhang 505582bec1SHong Zhang } PC_ML; 515582bec1SHong Zhang extern int PetscML_getrow(void *ML_data,int N_requested_rows,int requested_rows[], 525582bec1SHong Zhang int allocated_space,int columns[],double values[],int row_lengths[]); 535582bec1SHong Zhang extern int PetscML_matvec(void *ML_data, int in_length, double p[], int out_length,double ap[]); 545582bec1SHong Zhang extern int PetscML_comm(double x[], void *ML_data); 555582bec1SHong Zhang extern PetscErrorCode MatMult_ML(Mat,Vec,Vec); 565582bec1SHong Zhang extern PetscErrorCode MatMultAdd_ML(Mat,Vec,Vec,Vec); 575582bec1SHong Zhang extern PetscErrorCode MatConvert_MPIAIJ_ML(Mat,const MatType,Mat*); 585582bec1SHong Zhang extern PetscErrorCode MatDestroy_ML(Mat); 595582bec1SHong Zhang extern PetscErrorCode MatConvert_ML_SeqAIJ(FineGridCtx*,Mat*); 605582bec1SHong Zhang extern PetscErrorCode MatConvert_ML_SHELL(ML_Operator*,Mat*); 615582bec1SHong Zhang 625582bec1SHong Zhang /* -------------------------------------------------------------------------- */ 635582bec1SHong Zhang /* 645582bec1SHong Zhang PCSetUp_ML - Prepares for the use of the ML preconditioner 655582bec1SHong Zhang by setting data structures and options. 665582bec1SHong Zhang 675582bec1SHong Zhang Input Parameter: 685582bec1SHong Zhang . pc - the preconditioner context 695582bec1SHong Zhang 705582bec1SHong Zhang Application Interface Routine: PCSetUp() 715582bec1SHong Zhang 725582bec1SHong Zhang Notes: 735582bec1SHong Zhang The interface routine PCSetUp() is not usually called directly by 745582bec1SHong Zhang the user, but instead is called by PCApply() if necessary. 755582bec1SHong Zhang */ 765582bec1SHong Zhang 775582bec1SHong Zhang #undef __FUNCT__ 785582bec1SHong Zhang #define __FUNCT__ "PCSetUp_ML" 795582bec1SHong Zhang static PetscErrorCode PCSetUp_ML(PC pc) 805582bec1SHong Zhang { 815582bec1SHong Zhang PetscErrorCode ierr; 825582bec1SHong Zhang PetscMPIInt size,rank; 835582bec1SHong Zhang FineGridCtx *PetscMLdata; 845582bec1SHong Zhang ML *ml_object; 855582bec1SHong Zhang ML_Aggregate *agg_object; 865582bec1SHong Zhang ML_Operator *mlmat; 875582bec1SHong Zhang PetscInt nlocal_allcols,Nlevels,mllevel,level,m,fine_level; 885582bec1SHong Zhang Mat A,Aloc; 895582bec1SHong Zhang GridCtx *gridctx; 905582bec1SHong Zhang PC pc_coarse; 915582bec1SHong Zhang PC_ML *pc_ml=PETSC_NULL; 925582bec1SHong Zhang PetscObjectContainer container; 935582bec1SHong Zhang 945582bec1SHong Zhang PetscFunctionBegin; 955582bec1SHong Zhang ierr = PetscObjectQuery((PetscObject)pc,"PC_ML",(PetscObject *)&container);CHKERRQ(ierr); 965582bec1SHong Zhang if (container) { 975582bec1SHong Zhang ierr = PetscObjectContainerGetPointer(container,(void **)&pc_ml);CHKERRQ(ierr); 985582bec1SHong Zhang } else { 995582bec1SHong Zhang SETERRQ(PETSC_ERR_ARG_NULL,"Container does not exit"); 1005582bec1SHong Zhang } 1015582bec1SHong Zhang 1025582bec1SHong Zhang /* setup special features of PCML */ 1035582bec1SHong Zhang /*--------------------------------*/ 1045582bec1SHong Zhang /* covert A to Aloc to be used by ML at fine grid */ 1055582bec1SHong Zhang A = pc->pmat; 1065582bec1SHong Zhang ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr); 1075582bec1SHong Zhang ierr = MPI_Comm_rank(A->comm,&rank);CHKERRQ(ierr); 1085582bec1SHong Zhang pc_ml->size = size; 1095582bec1SHong Zhang if (size > 1){ 1105582bec1SHong Zhang Aloc = PETSC_NULL; 1115582bec1SHong Zhang ierr = MatConvert_MPIAIJ_ML(A,0,&Aloc);CHKERRQ(ierr); 1125582bec1SHong Zhang } else { 1135582bec1SHong Zhang Aloc = A; 1145582bec1SHong Zhang } 1155582bec1SHong Zhang 1165582bec1SHong Zhang /* create and initialize struct 'PetscMLdata' */ 1175582bec1SHong Zhang ierr = PetscNew(FineGridCtx,&PetscMLdata);CHKERRQ(ierr); 1185582bec1SHong Zhang ierr = PetscMemzero(PetscMLdata,sizeof(FineGridCtx));CHKERRQ(ierr); 1195582bec1SHong Zhang PetscMLdata->A = A; 1205582bec1SHong Zhang PetscMLdata->Aloc = Aloc; 1215582bec1SHong Zhang ierr = PetscMalloc((Aloc->n+1)*sizeof(PetscScalar),&PetscMLdata->pwork);CHKERRQ(ierr); 1225582bec1SHong Zhang pc_ml->PetscMLdata = PetscMLdata; 1235582bec1SHong Zhang 1245582bec1SHong Zhang /* create ML discretization matrix at fine grid */ 1255582bec1SHong Zhang ierr = MatGetSize(Aloc,&m,&nlocal_allcols);CHKERRQ(ierr); 1265582bec1SHong Zhang ML_Create(&ml_object,pc_ml->MaxNlevels); 1275582bec1SHong Zhang ML_Init_Amatrix(ml_object,0,m,m,PetscMLdata); 1285582bec1SHong Zhang ML_Set_Amatrix_Getrow(ml_object,0,PetscML_getrow,PetscML_comm,nlocal_allcols); 1295582bec1SHong Zhang ML_Set_Amatrix_Matvec(ml_object,0,PetscML_matvec); 1305582bec1SHong Zhang 1315582bec1SHong Zhang /* aggregation */ 1325582bec1SHong Zhang ML_Aggregate_Create(&agg_object); 1335582bec1SHong Zhang ML_Aggregate_Set_MaxCoarseSize(agg_object,pc_ml->MaxCoarseSize); 1345582bec1SHong Zhang /* set options */ 1355582bec1SHong Zhang switch (pc_ml->CoarsenScheme) { 1365582bec1SHong Zhang case 1: 1375582bec1SHong Zhang ML_Aggregate_Set_CoarsenScheme_Coupled(agg_object);break; 1385582bec1SHong Zhang case 2: 1395582bec1SHong Zhang ML_Aggregate_Set_CoarsenScheme_MIS(agg_object);break; 1405582bec1SHong Zhang case 3: 1415582bec1SHong Zhang ML_Aggregate_Set_CoarsenScheme_METIS(agg_object);break; 1425582bec1SHong Zhang } 1435582bec1SHong Zhang ML_Aggregate_Set_Threshold(agg_object,pc_ml->Threshold); 1445582bec1SHong Zhang ML_Aggregate_Set_DampingFactor(agg_object,pc_ml->DampingFactor); 1455582bec1SHong Zhang if (pc_ml->SpectralNormScheme_Anorm){ 1465582bec1SHong Zhang ML_Aggregate_Set_SpectralNormScheme_Anorm(agg_object); 1475582bec1SHong Zhang } 1485582bec1SHong Zhang 1495582bec1SHong Zhang Nlevels = ML_Gen_MGHierarchy_UsingAggregation(ml_object,0,ML_INCREASING,agg_object); 1505582bec1SHong Zhang if (Nlevels<=0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Nlevels %d must > 0",Nlevels); 1515582bec1SHong Zhang ierr = MGSetLevels(pc,Nlevels,PETSC_NULL);CHKERRQ(ierr); 1525582bec1SHong Zhang pc_ml->ml_object = ml_object; 1535582bec1SHong Zhang pc_ml->agg_object = agg_object; 1545582bec1SHong Zhang 1555582bec1SHong Zhang ierr = PetscMalloc(Nlevels*sizeof(GridCtx),&gridctx);CHKERRQ(ierr); 1565582bec1SHong Zhang fine_level = Nlevels - 1; 1575582bec1SHong Zhang pc_ml->gridctx = gridctx; 1585582bec1SHong Zhang pc_ml->fine_level = fine_level; 1595582bec1SHong Zhang 1605582bec1SHong Zhang /* wrap ML matrices by PETSc shell matrices at coarsened grids. 1615582bec1SHong Zhang Level 0 is the finest grid for ML, but coarsest for PETSc! */ 1625582bec1SHong Zhang gridctx[fine_level].A = A; 1635582bec1SHong Zhang level = fine_level - 1; 1645582bec1SHong Zhang if (size == 1){ /* convert ML mat format into petsc seqaij format */ 1655582bec1SHong Zhang PetscMLdata->rlen_max = A->N; 1665582bec1SHong Zhang ierr = PetscMalloc(PetscMLdata->rlen_max*(sizeof(PetscScalar)+sizeof(PetscInt)),&PetscMLdata->cols);CHKERRQ(ierr); 1675582bec1SHong Zhang PetscMLdata->vals = (PetscScalar*)(PetscMLdata->cols + PetscMLdata->rlen_max); 1685582bec1SHong Zhang for (mllevel=1; mllevel<Nlevels; mllevel++){ 1695582bec1SHong Zhang PetscMLdata->mlmat = &(ml_object->Pmat[mllevel]); 1705582bec1SHong Zhang ierr = MatConvert_ML_SeqAIJ(PetscMLdata,&gridctx[level].P);CHKERRQ(ierr); 1715582bec1SHong Zhang PetscMLdata->mlmat = &(ml_object->Amat[mllevel]); 1725582bec1SHong Zhang ierr = MatConvert_ML_SeqAIJ(PetscMLdata,&gridctx[level].A);CHKERRQ(ierr); 1735582bec1SHong Zhang PetscMLdata->mlmat = &(ml_object->Rmat[mllevel-1]); 1745582bec1SHong Zhang ierr = MatConvert_ML_SeqAIJ(PetscMLdata,&gridctx[level].R);CHKERRQ(ierr); 1755582bec1SHong Zhang level--; 1765582bec1SHong Zhang } 1775582bec1SHong Zhang } else { /* convert ML mat format into petsc shell format */ 1785582bec1SHong Zhang for (mllevel=1; mllevel<Nlevels; mllevel++){ 1795582bec1SHong Zhang mlmat = &(ml_object->Pmat[mllevel]); 1805582bec1SHong Zhang ierr = MatConvert_ML_SHELL(mlmat,&gridctx[level].P);CHKERRQ(ierr); 1815582bec1SHong Zhang mlmat = &(ml_object->Amat[mllevel]); 1825582bec1SHong Zhang ierr = MatConvert_ML_SHELL(mlmat,&gridctx[level].A);CHKERRQ(ierr); 1835582bec1SHong Zhang mlmat = &(ml_object->Rmat[mllevel-1]); 1845582bec1SHong Zhang ierr = MatConvert_ML_SHELL(mlmat,&gridctx[level].R);CHKERRQ(ierr); 1855582bec1SHong Zhang level--; 1865582bec1SHong Zhang } 1875582bec1SHong Zhang } 1885582bec1SHong Zhang 1895582bec1SHong Zhang /* create coarse level and the interpolation between the levels */ 1905582bec1SHong Zhang level = fine_level; 1915582bec1SHong Zhang while ( level >= 0 ){ 1925582bec1SHong Zhang if (level != fine_level){ 1935582bec1SHong Zhang ierr = VecCreate(gridctx[level].A->comm,&gridctx[level].x);CHKERRQ(ierr); 1945582bec1SHong Zhang ierr = VecSetSizes(gridctx[level].x,gridctx[level].A->n,PETSC_DECIDE);CHKERRQ(ierr); 1955582bec1SHong Zhang ierr = VecSetType(gridctx[level].x,VECMPI);CHKERRQ(ierr); 1965582bec1SHong Zhang ierr = MGSetX(pc,level,gridctx[level].x);CHKERRQ(ierr); 1975582bec1SHong Zhang 1985582bec1SHong Zhang ierr = VecCreate(gridctx[level].A->comm,&gridctx[level].b);CHKERRQ(ierr); 1995582bec1SHong Zhang ierr = VecSetSizes(gridctx[level].b,gridctx[level].A->m,PETSC_DECIDE);CHKERRQ(ierr); 2005582bec1SHong Zhang ierr = VecSetType(gridctx[level].b,VECMPI);CHKERRQ(ierr); 2015582bec1SHong Zhang ierr = MGSetRhs(pc,level,gridctx[level].b);CHKERRQ(ierr); 2025582bec1SHong Zhang } 2035582bec1SHong Zhang ierr = VecCreate(gridctx[level].A->comm,&gridctx[level].r);CHKERRQ(ierr); 2045582bec1SHong Zhang ierr = VecSetSizes(gridctx[level].r,gridctx[level].A->m,PETSC_DECIDE);CHKERRQ(ierr); 2055582bec1SHong Zhang ierr = VecSetType(gridctx[level].r,VECMPI);CHKERRQ(ierr); 2065582bec1SHong Zhang ierr = MGSetR(pc,level,gridctx[level].r);CHKERRQ(ierr); 2075582bec1SHong Zhang 2085582bec1SHong Zhang if (level == 0){ 2095582bec1SHong Zhang ierr = MGGetCoarseSolve(pc,&gridctx[level].ksp);CHKERRQ(ierr); 2105582bec1SHong Zhang } else { 2115582bec1SHong Zhang ierr = MGGetSmoother(pc,level,&gridctx[level].ksp);CHKERRQ(ierr); 2125582bec1SHong Zhang ierr = MGSetResidual(pc,level,MGDefaultResidual,gridctx[level].A);CHKERRQ(ierr); 2135582bec1SHong Zhang if (level == fine_level){ 2145582bec1SHong Zhang ierr = KSPSetOptionsPrefix(gridctx[level].ksp,"mg_fine_");CHKERRQ(ierr); 2155582bec1SHong Zhang ierr = MGSetR(pc,level,gridctx[level].r);CHKERRQ(ierr); 2165582bec1SHong Zhang } 2175582bec1SHong Zhang } 2185582bec1SHong Zhang ierr = KSPSetOperators(gridctx[level].ksp,gridctx[level].A,gridctx[level].A,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 2195582bec1SHong Zhang 2205582bec1SHong Zhang if (level < fine_level){ 2215582bec1SHong Zhang if (size > 1){ 2225582bec1SHong Zhang ierr = KSPGetPC(gridctx[level].ksp,&pc_coarse);CHKERRQ(ierr); 2235582bec1SHong Zhang ierr = PCSetType(pc_coarse,PCNONE);CHKERRQ(ierr); 2245582bec1SHong Zhang } 2255582bec1SHong Zhang ierr = MGSetInterpolate(pc,level+1,gridctx[level].P);CHKERRQ(ierr); 2265582bec1SHong Zhang ierr = MGSetRestriction(pc,level+1,gridctx[level].R);CHKERRQ(ierr); 2275582bec1SHong Zhang } 2285582bec1SHong Zhang level--; 2295582bec1SHong Zhang } 2305582bec1SHong Zhang 2315582bec1SHong Zhang /* now call PCSetUp_MG() */ 2325582bec1SHong Zhang /*--------------------------------*/ 2335582bec1SHong Zhang ierr = (*pc_ml->PCSetUp)(pc);CHKERRQ(ierr); 2345582bec1SHong Zhang PetscFunctionReturn(0); 2355582bec1SHong Zhang } 2365582bec1SHong Zhang 2375582bec1SHong Zhang #undef __FUNCT__ 2385582bec1SHong Zhang #define __FUNCT__ "PetscObjectContainerDestroy_PC_ML" 2395582bec1SHong Zhang PetscErrorCode PetscObjectContainerDestroy_PC_ML(void *ptr) 2405582bec1SHong Zhang { 2415582bec1SHong Zhang PetscErrorCode ierr; 2425582bec1SHong Zhang PC_ML *pc_ml = (PC_ML*)ptr; 2435582bec1SHong Zhang PetscInt level; 2445582bec1SHong Zhang 2455582bec1SHong Zhang PetscFunctionBegin; 2465582bec1SHong Zhang if (pc_ml->size > 1){ierr = MatDestroy(pc_ml->PetscMLdata->Aloc);CHKERRQ(ierr);} 2475582bec1SHong Zhang ML_Aggregate_Destroy(&pc_ml->agg_object); 2485582bec1SHong Zhang ML_Destroy(&pc_ml->ml_object); 2495582bec1SHong Zhang 2505582bec1SHong Zhang ierr = PetscFree(pc_ml->PetscMLdata->pwork);CHKERRQ(ierr); 2515582bec1SHong Zhang if (pc_ml->size == 1){ierr = PetscFree(pc_ml->PetscMLdata->cols);CHKERRQ(ierr);} 2525582bec1SHong Zhang ierr = PetscFree(pc_ml->PetscMLdata);CHKERRQ(ierr); 2535582bec1SHong Zhang 2545582bec1SHong Zhang level = pc_ml->fine_level; 2555582bec1SHong Zhang while ( level>= 0 ){ 2565582bec1SHong Zhang if (level != pc_ml->fine_level){ 2575582bec1SHong Zhang ierr = MatDestroy(pc_ml->gridctx[level].A);CHKERRQ(ierr); 2585582bec1SHong Zhang ierr = MatDestroy(pc_ml->gridctx[level].P);CHKERRQ(ierr); 2595582bec1SHong Zhang ierr = MatDestroy(pc_ml->gridctx[level].R);CHKERRQ(ierr); 2605582bec1SHong Zhang ierr = VecDestroy(pc_ml->gridctx[level].x);CHKERRQ(ierr); 2615582bec1SHong Zhang ierr = VecDestroy(pc_ml->gridctx[level].b);CHKERRQ(ierr); 2625582bec1SHong Zhang } 2635582bec1SHong Zhang ierr = VecDestroy(pc_ml->gridctx[level].r);CHKERRQ(ierr); 2645582bec1SHong Zhang level--; 2655582bec1SHong Zhang } 2665582bec1SHong Zhang ierr = PetscFree(pc_ml->gridctx);CHKERRQ(ierr); 2675582bec1SHong Zhang ierr = PetscFree(pc_ml);CHKERRQ(ierr); 2685582bec1SHong Zhang PetscFunctionReturn(0); 2695582bec1SHong Zhang } 2705582bec1SHong Zhang /* -------------------------------------------------------------------------- */ 2715582bec1SHong Zhang /* 2725582bec1SHong Zhang PCDestroy_ML - Destroys the private context for the ML preconditioner 2735582bec1SHong Zhang that was created with PCCreate_ML(). 2745582bec1SHong Zhang 2755582bec1SHong Zhang Input Parameter: 2765582bec1SHong Zhang . pc - the preconditioner context 2775582bec1SHong Zhang 2785582bec1SHong Zhang Application Interface Routine: PCDestroy() 2795582bec1SHong Zhang */ 2805582bec1SHong Zhang #undef __FUNCT__ 2815582bec1SHong Zhang #define __FUNCT__ "PCDestroy_ML" 2825582bec1SHong Zhang static PetscErrorCode PCDestroy_ML(PC pc) 2835582bec1SHong Zhang { 2845582bec1SHong Zhang PetscErrorCode ierr; 2855582bec1SHong Zhang PC_ML *pc_ml=PETSC_NULL; 2865582bec1SHong Zhang PetscObjectContainer container; 2875582bec1SHong Zhang 2885582bec1SHong Zhang PetscFunctionBegin; 2895582bec1SHong Zhang ierr = PetscObjectQuery((PetscObject)pc,"PC_ML",(PetscObject *)&container);CHKERRQ(ierr); 2905582bec1SHong Zhang if (container) { 2915582bec1SHong Zhang ierr = PetscObjectContainerGetPointer(container,(void **)&pc_ml);CHKERRQ(ierr); 2925582bec1SHong Zhang pc->ops->destroy = pc_ml->PCDestroy; 2935582bec1SHong Zhang } else { 2945582bec1SHong Zhang SETERRQ(PETSC_ERR_ARG_NULL,"Container does not exit"); 2955582bec1SHong Zhang } 296*9cb0ec93SHong Zhang /* detach pc and PC_ML and dereference container */ 2975582bec1SHong Zhang ierr = PetscObjectCompose((PetscObject)pc,"PC_ML",0);CHKERRQ(ierr); 2985582bec1SHong Zhang ierr = (*pc->ops->destroy)(pc);CHKERRQ(ierr); 2995582bec1SHong Zhang 3005582bec1SHong Zhang ierr = PetscObjectContainerDestroy(container);CHKERRQ(ierr); 3015582bec1SHong Zhang PetscFunctionReturn(0); 3025582bec1SHong Zhang } 3035582bec1SHong Zhang 3045582bec1SHong Zhang #undef __FUNCT__ 3055582bec1SHong Zhang #define __FUNCT__ "PCSetFromOptions_ML" 3065582bec1SHong Zhang static PetscErrorCode PCSetFromOptions_ML(PC pc) 3075582bec1SHong Zhang { 3085582bec1SHong Zhang PetscErrorCode ierr; 3095582bec1SHong Zhang PetscInt indx,m,PrintLevel,MaxNlevels,MaxCoarseSize; 3105582bec1SHong Zhang PetscReal Threshold,DampingFactor; 3115582bec1SHong Zhang PetscTruth flg; 3125582bec1SHong Zhang const char *type[] = {"additive","multiplicative","full","cascade","kascade"}; 3135582bec1SHong Zhang const char *scheme[] = {"Uncoupled","Coupled","MIS","METIS"}; 3145582bec1SHong Zhang PC_ML *pc_ml=PETSC_NULL; 3155582bec1SHong Zhang PetscObjectContainer container; 3165582bec1SHong Zhang 3175582bec1SHong Zhang PetscFunctionBegin; 3185582bec1SHong Zhang ierr = PetscObjectQuery((PetscObject)pc,"PC_ML",(PetscObject *)&container);CHKERRQ(ierr); 3195582bec1SHong Zhang if (container) { 3205582bec1SHong Zhang ierr = PetscObjectContainerGetPointer(container,(void **)&pc_ml);CHKERRQ(ierr); 3215582bec1SHong Zhang } else { 3225582bec1SHong Zhang SETERRQ(PETSC_ERR_ARG_NULL,"Container does not exit"); 3235582bec1SHong Zhang } 3245582bec1SHong Zhang ierr = PetscOptionsHead("MG options");CHKERRQ(ierr); 3255582bec1SHong Zhang /* inherited MG options */ 3265582bec1SHong Zhang ierr = PetscOptionsInt("-pc_mg_cycles","1 for V cycle, 2 for W-cycle","MGSetCycles",1,&m,&flg);CHKERRQ(ierr); 3275582bec1SHong Zhang if (flg) { 3285582bec1SHong Zhang ierr = MGSetCycles(pc,m);CHKERRQ(ierr); 3295582bec1SHong Zhang } 3305582bec1SHong Zhang ierr = PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","MGSetNumberSmoothUp",1,&m,&flg);CHKERRQ(ierr); 3315582bec1SHong Zhang if (flg) { 3325582bec1SHong Zhang ierr = MGSetNumberSmoothUp(pc,m);CHKERRQ(ierr); 3335582bec1SHong Zhang } 3345582bec1SHong Zhang ierr = PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","MGSetNumberSmoothDown",1,&m,&flg);CHKERRQ(ierr); 3355582bec1SHong Zhang if (flg) { 3365582bec1SHong Zhang ierr = MGSetNumberSmoothDown(pc,m);CHKERRQ(ierr); 3375582bec1SHong Zhang } 3385582bec1SHong Zhang ierr = PetscOptionsEList("-pc_mg_type","Multigrid type","MGSetType",type,5,type[1],&indx,&flg);CHKERRQ(ierr); 3395582bec1SHong Zhang if (flg) { 3405582bec1SHong Zhang MGType mg = (MGType) 0; 3415582bec1SHong Zhang switch (indx) { 3425582bec1SHong Zhang case 0: 3435582bec1SHong Zhang mg = MGADDITIVE; 3445582bec1SHong Zhang break; 3455582bec1SHong Zhang case 1: 3465582bec1SHong Zhang mg = MGMULTIPLICATIVE; 3475582bec1SHong Zhang break; 3485582bec1SHong Zhang case 2: 3495582bec1SHong Zhang mg = MGFULL; 3505582bec1SHong Zhang break; 3515582bec1SHong Zhang case 3: 3525582bec1SHong Zhang mg = MGKASKADE; 3535582bec1SHong Zhang break; 3545582bec1SHong Zhang case 4: 3555582bec1SHong Zhang mg = MGKASKADE; 3565582bec1SHong Zhang break; 3575582bec1SHong Zhang } 3585582bec1SHong Zhang ierr = MGSetType(pc,mg);CHKERRQ(ierr); 3595582bec1SHong Zhang } 3605582bec1SHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 3615582bec1SHong Zhang 3625582bec1SHong Zhang /* ML options */ 3635582bec1SHong Zhang ierr = PetscOptionsHead("ML options");CHKERRQ(ierr); 3645582bec1SHong Zhang /* set defaults */ 3655582bec1SHong Zhang PrintLevel = 0; 3665582bec1SHong Zhang MaxNlevels = 10; 3675582bec1SHong Zhang MaxCoarseSize = 1; 3685582bec1SHong Zhang indx = 0; 3695582bec1SHong Zhang Threshold = 0.0; 3705582bec1SHong Zhang DampingFactor = 4.0/3.0; 3715582bec1SHong Zhang 3725582bec1SHong Zhang ierr = PetscOptionsInt("-pc_ml_PrintLevel","Print level","ML_Set_PrintLevel",PrintLevel,&PrintLevel,PETSC_NULL);CHKERRQ(ierr); 3735582bec1SHong Zhang ML_Set_PrintLevel(PrintLevel); 3745582bec1SHong Zhang 3755582bec1SHong Zhang ierr = PetscOptionsInt("-pc_ml_maxNlevels","Maximum number of levels","None",MaxNlevels,&MaxNlevels,PETSC_NULL);CHKERRQ(ierr); 3765582bec1SHong Zhang pc_ml->MaxNlevels = MaxNlevels; 3775582bec1SHong Zhang 3785582bec1SHong Zhang ierr = PetscOptionsInt("-pc_ml_maxCoarseSize","Maximum coarsest mesh size","ML_Aggregate_Set_MaxCoarseSize",MaxCoarseSize,&MaxCoarseSize,PETSC_NULL);CHKERRQ(ierr); 3795582bec1SHong Zhang pc_ml->MaxCoarseSize = MaxCoarseSize; 3805582bec1SHong Zhang 3815582bec1SHong Zhang ierr = PetscOptionsEList("-pc_ml_CoarsenScheme","Aggregate Coarsen Scheme","ML_Aggregate_Set_CoarsenScheme_*",scheme,4,scheme[0],&indx,PETSC_NULL);CHKERRQ(ierr); 3825582bec1SHong Zhang pc_ml->CoarsenScheme = indx; 3835582bec1SHong Zhang 3845582bec1SHong Zhang ierr = PetscOptionsReal("-pc_ml_DampingFactor","P damping factor","ML_Aggregate_Set_DampingFactor",DampingFactor,&DampingFactor,PETSC_NULL);CHKERRQ(ierr); 3855582bec1SHong Zhang pc_ml->DampingFactor = DampingFactor; 3865582bec1SHong Zhang 3875582bec1SHong Zhang ierr = PetscOptionsReal("-pc_ml_Threshold","Smoother drop tol","ML_Aggregate_Set_Threshold",Threshold,&Threshold,PETSC_NULL);CHKERRQ(ierr); 3885582bec1SHong Zhang pc_ml->Threshold = Threshold; 3895582bec1SHong Zhang 3905582bec1SHong Zhang ierr = PetscOptionsLogical("-pc_ml_SpectralNormScheme_Anorm","Method used for estimating spectral radius","ML_Aggregate_Set_SpectralNormScheme_Anorm",PETSC_FALSE,&pc_ml->SpectralNormScheme_Anorm,PETSC_FALSE); 3915582bec1SHong Zhang 3925582bec1SHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 3935582bec1SHong Zhang PetscFunctionReturn(0); 3945582bec1SHong Zhang } 3955582bec1SHong Zhang 3965582bec1SHong Zhang /* -------------------------------------------------------------------------- */ 3975582bec1SHong Zhang /* 3985582bec1SHong Zhang PCCreate_ML - Creates a ML preconditioner context, PC_ML, 3995582bec1SHong Zhang and sets this as the private data within the generic preconditioning 4005582bec1SHong Zhang context, PC, that was created within PCCreate(). 4015582bec1SHong Zhang 4025582bec1SHong Zhang Input Parameter: 4035582bec1SHong Zhang . pc - the preconditioner context 4045582bec1SHong Zhang 4055582bec1SHong Zhang Application Interface Routine: PCCreate() 4065582bec1SHong Zhang */ 4075582bec1SHong Zhang 4085582bec1SHong Zhang /*MC 4095582bec1SHong Zhang PCML - Use geometric multigrid preconditioning. This preconditioner requires you provide 4105582bec1SHong Zhang fine grid discretization matrix. The coarser grid matrices and restriction/interpolation 4115582bec1SHong Zhang operators are computed by ML and wrapped as PETSc shell matrices. 4125582bec1SHong Zhang 4135582bec1SHong Zhang Options Database Key: (not done yet!) 4145582bec1SHong Zhang + -pc_mg_maxlevels <nlevels> - maximum number of levels including finest 4155582bec1SHong Zhang . -pc_mg_cycles 1 or 2 - for V or W-cycle 4165582bec1SHong Zhang . -pc_mg_smoothup <n> - number of smoothing steps after interpolation 4175582bec1SHong Zhang . -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator 4185582bec1SHong Zhang . -pc_mg_type <additive,multiplicative,full,cascade> - multiplicative is the default 4195582bec1SHong Zhang . -pc_mg_monitor - print information on the multigrid convergence 4205582bec1SHong Zhang - -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices 4215582bec1SHong Zhang to the Socket viewer for reading from Matlab. 4225582bec1SHong Zhang 4235582bec1SHong Zhang Level: intermediate 4245582bec1SHong Zhang 4255582bec1SHong Zhang Concepts: multigrid 4265582bec1SHong Zhang 4275582bec1SHong Zhang .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, 4285582bec1SHong Zhang MGSetLevels(), MGGetLevels(), MGSetType(), MPSetCycles(), MGSetNumberSmoothDown(), 4295582bec1SHong Zhang MGSetNumberSmoothUp(), MGGetCoarseSolve(), MGSetResidual(), MGSetInterpolation(), 4305582bec1SHong Zhang MGSetRestriction(), MGGetSmoother(), MGGetSmootherUp(), MGGetSmootherDown(), 4315582bec1SHong Zhang MGSetCyclesOnLevel(), MGSetRhs(), MGSetX(), MGSetR() 4325582bec1SHong Zhang M*/ 4335582bec1SHong Zhang 4345582bec1SHong Zhang EXTERN_C_BEGIN 4355582bec1SHong Zhang #undef __FUNCT__ 4365582bec1SHong Zhang #define __FUNCT__ "PCCreate_ML" 4375582bec1SHong Zhang PetscErrorCode PCCreate_ML(PC pc) 4385582bec1SHong Zhang { 4395582bec1SHong Zhang PetscErrorCode ierr; 4405582bec1SHong Zhang PC_ML *pc_ml; 4415582bec1SHong Zhang PetscObjectContainer container; 4425582bec1SHong Zhang 4435582bec1SHong Zhang PetscFunctionBegin; 4445582bec1SHong Zhang /* initialize pc as PCMG */ 4455582bec1SHong Zhang ierr = PCSetType(pc,PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */ 4465582bec1SHong Zhang 4475582bec1SHong Zhang /* create a supporting struct and attach it to pc */ 4485582bec1SHong Zhang ierr = PetscNew(PC_ML,&pc_ml);CHKERRQ(ierr); 449*9cb0ec93SHong Zhang ierr = PetscMemzero(pc_ml,sizeof(PC_ML));CHKERRQ(ierr); 4505582bec1SHong Zhang ierr = PetscObjectContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr); 4515582bec1SHong Zhang ierr = PetscObjectContainerSetPointer(container,pc_ml);CHKERRQ(ierr); 452*9cb0ec93SHong Zhang ierr = PetscObjectContainerSetUserDestroy(container,PetscObjectContainerDestroy_PC_ML);CHKERRQ(ierr); 4535582bec1SHong Zhang ierr = PetscObjectCompose((PetscObject)pc,"PC_ML",(PetscObject)container);CHKERRQ(ierr); 4545582bec1SHong Zhang 4555582bec1SHong Zhang pc_ml->PCSetUp = pc->ops->setup; 4565582bec1SHong Zhang pc_ml->PCDestroy = pc->ops->destroy; 4575582bec1SHong Zhang 4585582bec1SHong Zhang /* overwrite the pointers of PCMG by the functions of PCML */ 4595582bec1SHong Zhang pc->ops->setfromoptions = PCSetFromOptions_ML; 4605582bec1SHong Zhang pc->ops->setup = PCSetUp_ML; 4615582bec1SHong Zhang pc->ops->destroy = PCDestroy_ML; 4625582bec1SHong Zhang PetscFunctionReturn(0); 4635582bec1SHong Zhang } 4645582bec1SHong Zhang EXTERN_C_END 4655582bec1SHong Zhang 4665582bec1SHong Zhang int PetscML_getrow(void *ML_data, int N_requested_rows, int requested_rows[], 4675582bec1SHong Zhang int allocated_space, int columns[], double values[], int row_lengths[]) 4685582bec1SHong Zhang { 4695582bec1SHong Zhang PetscErrorCode ierr; 4705582bec1SHong Zhang Mat Aloc; 4715582bec1SHong Zhang Mat_SeqAIJ *a; 4725582bec1SHong Zhang PetscInt m,i,j,k=0,row,*aj; 4735582bec1SHong Zhang PetscScalar *aa; 4745582bec1SHong Zhang FineGridCtx *ml=(FineGridCtx*)ML_data; 4755582bec1SHong Zhang 4765582bec1SHong Zhang Aloc = ml->Aloc; 4775582bec1SHong Zhang a = (Mat_SeqAIJ*)Aloc->data; 4785582bec1SHong Zhang ierr = MatGetSize(Aloc,&m,PETSC_NULL);CHKERRQ(ierr); 4795582bec1SHong Zhang 4805582bec1SHong Zhang for (i = 0; i<N_requested_rows; i++) { 4815582bec1SHong Zhang row = requested_rows[i]; 4825582bec1SHong Zhang row_lengths[i] = a->ilen[row]; 4835582bec1SHong Zhang if (allocated_space < k+row_lengths[i]) return(0); 4845582bec1SHong Zhang if ( (row >= 0) || (row <= (m-1)) ) { 4855582bec1SHong Zhang aj = a->j + a->i[row]; 4865582bec1SHong Zhang aa = a->a + a->i[row]; 4875582bec1SHong Zhang for (j=0; j<row_lengths[i]; j++){ 4885582bec1SHong Zhang columns[k] = aj[j]; 4895582bec1SHong Zhang values[k++] = aa[j]; 4905582bec1SHong Zhang } 4915582bec1SHong Zhang } 4925582bec1SHong Zhang } 4935582bec1SHong Zhang return(1); 4945582bec1SHong Zhang } 4955582bec1SHong Zhang 4965582bec1SHong Zhang int PetscML_matvec(void *ML_data,int in_length,double p[],int out_length,double ap[]) 4975582bec1SHong Zhang { 4985582bec1SHong Zhang PetscErrorCode ierr; 4995582bec1SHong Zhang FineGridCtx *ml=(FineGridCtx*)ML_data; 5005582bec1SHong Zhang Mat A=ml->A, Aloc=ml->Aloc; 5015582bec1SHong Zhang Vec x,y; 5025582bec1SHong Zhang PetscMPIInt size; 5035582bec1SHong Zhang PetscScalar *pwork=ml->pwork; 5045582bec1SHong Zhang PetscInt i; 5055582bec1SHong Zhang 5065582bec1SHong Zhang ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr); 5075582bec1SHong Zhang ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,out_length,ap,&y);CHKERRQ(ierr); 5085582bec1SHong Zhang if (size == 1){ 5095582bec1SHong Zhang ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,in_length,p,&x);CHKERRQ(ierr); 5105582bec1SHong Zhang } else { 5115582bec1SHong Zhang for (i=0; i<in_length; i++) pwork[i] = p[i]; 5125582bec1SHong Zhang PetscML_comm(pwork,ml); 5135582bec1SHong Zhang ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,Aloc->n,pwork,&x);CHKERRQ(ierr); 5145582bec1SHong Zhang } 5155582bec1SHong Zhang ierr = MatMult(Aloc,x,y);CHKERRQ(ierr); 5165582bec1SHong Zhang ierr = VecDestroy(x);CHKERRQ(ierr); 5175582bec1SHong Zhang ierr = VecDestroy(y);CHKERRQ(ierr); 5185582bec1SHong Zhang return 0; 5195582bec1SHong Zhang } 5205582bec1SHong Zhang 5215582bec1SHong Zhang int PetscML_comm(double p[],void *ML_data) 5225582bec1SHong Zhang { 5235582bec1SHong Zhang PetscErrorCode ierr; 5245582bec1SHong Zhang FineGridCtx *ml=(FineGridCtx*)ML_data; 5255582bec1SHong Zhang Mat A=ml->A; 5265582bec1SHong Zhang Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data; 5275582bec1SHong Zhang Vec x; 5285582bec1SHong Zhang PetscMPIInt size; 5295582bec1SHong Zhang PetscInt i,in_length=A->m,out_length=ml->Aloc->n; 5305582bec1SHong Zhang PetscScalar *array; 5315582bec1SHong Zhang 5325582bec1SHong Zhang ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr); 5335582bec1SHong Zhang if (size == 1) return 0; 5345582bec1SHong Zhang ierr = VecCreateSeqWithArray(PETSC_COMM_SELF,in_length,p,&x);CHKERRQ(ierr); 5355582bec1SHong Zhang ierr = VecScatterBegin(x,a->lvec,INSERT_VALUES,SCATTER_FORWARD,a->Mvctx);CHKERRQ(ierr); 5365582bec1SHong Zhang ierr = VecScatterEnd(x,a->lvec,INSERT_VALUES,SCATTER_FORWARD,a->Mvctx);CHKERRQ(ierr); 5375582bec1SHong Zhang ierr = VecGetArray(a->lvec,&array);CHKERRQ(ierr); 5385582bec1SHong Zhang for (i=in_length; i<out_length; i++){ 5395582bec1SHong Zhang p[i] = array[i-in_length]; 5405582bec1SHong Zhang } 5415582bec1SHong Zhang ierr = VecDestroy(x);CHKERRQ(ierr); 5425582bec1SHong Zhang return 0; 5435582bec1SHong Zhang } 5445582bec1SHong Zhang #undef __FUNCT__ 5455582bec1SHong Zhang #define __FUNCT__ "MatMult_ML" 5465582bec1SHong Zhang PetscErrorCode MatMult_ML(Mat A,Vec x,Vec y) 5475582bec1SHong Zhang { 5485582bec1SHong Zhang PetscErrorCode ierr; 5495582bec1SHong Zhang Mat_MLShell *shell; 5505582bec1SHong Zhang PetscScalar *xarray,*yarray; 5515582bec1SHong Zhang PetscInt x_length,y_length; 5525582bec1SHong Zhang 5535582bec1SHong Zhang PetscFunctionBegin; 5545582bec1SHong Zhang ierr = MatShellGetContext(A,(void *)&shell);CHKERRQ(ierr); 5555582bec1SHong Zhang ierr = VecGetArray(x,&xarray);CHKERRQ(ierr); 5565582bec1SHong Zhang ierr = VecGetArray(y,&yarray);CHKERRQ(ierr); 5575582bec1SHong Zhang x_length = shell->mlmat->invec_leng; 5585582bec1SHong Zhang y_length = shell->mlmat->outvec_leng; 5595582bec1SHong Zhang 5605582bec1SHong Zhang ML_Operator_Apply(shell->mlmat,x_length,xarray,y_length,yarray); 5615582bec1SHong Zhang 5625582bec1SHong Zhang ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr); 5635582bec1SHong Zhang ierr = VecRestoreArray(y,&yarray);CHKERRQ(ierr); 5645582bec1SHong Zhang PetscFunctionReturn(0); 5655582bec1SHong Zhang } 5665582bec1SHong Zhang /* MatMultAdd_ML - Compute y = w + A*x */ 5675582bec1SHong Zhang #undef __FUNCT__ 5685582bec1SHong Zhang #define __FUNCT__ "MatMultAdd_ML" 5695582bec1SHong Zhang PetscErrorCode MatMultAdd_ML(Mat A,Vec x,Vec w,Vec y) 5705582bec1SHong Zhang { 5715582bec1SHong Zhang PetscErrorCode ierr; 5725582bec1SHong Zhang Mat_MLShell *shell; 5735582bec1SHong Zhang PetscScalar *xarray,*yarray; 5745582bec1SHong Zhang const PetscScalar one=1.0; 5755582bec1SHong Zhang PetscInt x_length,y_length; 5765582bec1SHong Zhang 5775582bec1SHong Zhang PetscFunctionBegin; 5785582bec1SHong Zhang ierr = MatShellGetContext(A,(void *)&shell);CHKERRQ(ierr); 5795582bec1SHong Zhang ierr = VecGetArray(x,&xarray);CHKERRQ(ierr); 5805582bec1SHong Zhang ierr = VecGetArray(y,&yarray);CHKERRQ(ierr); 5815582bec1SHong Zhang 5825582bec1SHong Zhang x_length = shell->mlmat->invec_leng; 5835582bec1SHong Zhang y_length = shell->mlmat->outvec_leng; 5845582bec1SHong Zhang 5855582bec1SHong Zhang ML_Operator_Apply(shell->mlmat,x_length,xarray,y_length,yarray); 5865582bec1SHong Zhang 5875582bec1SHong Zhang ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr); 5885582bec1SHong Zhang ierr = VecRestoreArray(y,&yarray);CHKERRQ(ierr); 5895582bec1SHong Zhang ierr = VecAXPY(&one,w,y);CHKERRQ(ierr); 5905582bec1SHong Zhang 5915582bec1SHong Zhang PetscFunctionReturn(0); 5925582bec1SHong Zhang } 5935582bec1SHong Zhang 5945582bec1SHong Zhang /* newtype is ignored because "ml" is not listed under Petsc MatType yet */ 5955582bec1SHong Zhang #undef __FUNCT__ 5965582bec1SHong Zhang #define __FUNCT__ "MatConvert_MPIAIJ_ML" 5975582bec1SHong Zhang PetscErrorCode MatConvert_MPIAIJ_ML(Mat A,const MatType newtype,Mat *Aloc) 5985582bec1SHong Zhang { 5995582bec1SHong Zhang PetscErrorCode ierr; 6005582bec1SHong Zhang Mat_MPIAIJ *mpimat=(Mat_MPIAIJ*)A->data; 6015582bec1SHong Zhang Mat_SeqAIJ *mat,*a=(Mat_SeqAIJ*)(mpimat->A)->data,*b=(Mat_SeqAIJ*)(mpimat->B)->data; 6025582bec1SHong Zhang PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j; 6035582bec1SHong Zhang PetscScalar *aa=a->a,*ba=b->a,*ca; 6045582bec1SHong Zhang PetscInt am=A->m,an=A->n,i,j,k; 6055582bec1SHong Zhang PetscInt *ci,*cj,ncols; 6065582bec1SHong Zhang MatReuse scall=MAT_INITIAL_MATRIX; 6075582bec1SHong Zhang 6085582bec1SHong Zhang PetscFunctionBegin; 6095582bec1SHong Zhang if (am != an) SETERRQ2(PETSC_ERR_ARG_WRONG,"A must have a square diagonal portion, am: %d != an: %d",am,an); 6105582bec1SHong Zhang 6115582bec1SHong Zhang if (*Aloc) scall = MAT_REUSE_MATRIX; 6125582bec1SHong Zhang if (scall == MAT_INITIAL_MATRIX){ 6135582bec1SHong Zhang ierr = PetscMalloc((1+am)*sizeof(PetscInt),&ci);CHKERRQ(ierr); 6145582bec1SHong Zhang ci[0] = 0; 6155582bec1SHong Zhang for (i=0; i<am; i++){ 6165582bec1SHong Zhang ci[i+1] = ci[i] + (ai[i+1] - ai[i]) + (bi[i+1] - bi[i]); 6175582bec1SHong Zhang } 6185582bec1SHong Zhang ierr = PetscMalloc((1+ci[am])*sizeof(PetscInt),&cj);CHKERRQ(ierr); 6195582bec1SHong Zhang ierr = PetscMalloc((1+ci[am])*sizeof(PetscScalar),&ca);CHKERRQ(ierr); 6205582bec1SHong Zhang 6215582bec1SHong Zhang k = 0; 6225582bec1SHong Zhang for (i=0; i<am; i++){ 6235582bec1SHong Zhang /* diagonal portion of A */ 6245582bec1SHong Zhang ncols = ai[i+1] - ai[i]; 6255582bec1SHong Zhang for (j=0; j<ncols; j++) { 6265582bec1SHong Zhang cj[k] = *aj++; 6275582bec1SHong Zhang ca[k++] = *aa++; 6285582bec1SHong Zhang } 6295582bec1SHong Zhang /* off-diagonal portion of A */ 6305582bec1SHong Zhang ncols = bi[i+1] - bi[i]; 6315582bec1SHong Zhang for (j=0; j<ncols; j++) { 6325582bec1SHong Zhang cj[k] = an + (*bj); bj++; 6335582bec1SHong Zhang ca[k++] = *ba++; 6345582bec1SHong Zhang } 6355582bec1SHong Zhang } 6365582bec1SHong Zhang if (k != ci[am]) SETERRQ2(PETSC_ERR_ARG_WRONG,"k: %d != ci[am]: %d",k,ci[am]); 6375582bec1SHong Zhang 6385582bec1SHong Zhang /* put together the new matrix */ 6395582bec1SHong Zhang an = mpimat->A->n+mpimat->B->n; 6405582bec1SHong Zhang ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,an,ci,cj,ca,Aloc);CHKERRQ(ierr); 6415582bec1SHong Zhang 6425582bec1SHong Zhang /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 6435582bec1SHong Zhang /* Since these are PETSc arrays, change flags to free them as necessary. */ 6445582bec1SHong Zhang mat = (Mat_SeqAIJ*)(*Aloc)->data; 6455582bec1SHong Zhang mat->freedata = PETSC_TRUE; 6465582bec1SHong Zhang mat->nonew = 0; 6475582bec1SHong Zhang } else if (scall == MAT_REUSE_MATRIX){ 6485582bec1SHong Zhang mat=(Mat_SeqAIJ*)(*Aloc)->data; 6495582bec1SHong Zhang ci = mat->i; cj = mat->j; ca = mat->a; 6505582bec1SHong Zhang for (i=0; i<am; i++) { 6515582bec1SHong Zhang /* diagonal portion of A */ 6525582bec1SHong Zhang ncols = ai[i+1] - ai[i]; 6535582bec1SHong Zhang for (j=0; j<ncols; j++) *ca++ = *aa++; 6545582bec1SHong Zhang /* off-diagonal portion of A */ 6555582bec1SHong Zhang ncols = bi[i+1] - bi[i]; 6565582bec1SHong Zhang for (j=0; j<ncols; j++) *ca++ = *ba++; 6575582bec1SHong Zhang } 6585582bec1SHong Zhang } else { 6595582bec1SHong Zhang SETERRQ1(PETSC_ERR_ARG_WRONG,"Invalid MatReuse %d",(int)scall); 6605582bec1SHong Zhang } 6615582bec1SHong Zhang PetscFunctionReturn(0); 6625582bec1SHong Zhang } 6635582bec1SHong Zhang extern PetscErrorCode MatDestroy_Shell(Mat); 6645582bec1SHong Zhang #undef __FUNCT__ 6655582bec1SHong Zhang #define __FUNCT__ "MatDestroy_ML" 6665582bec1SHong Zhang PetscErrorCode MatDestroy_ML(Mat A) 6675582bec1SHong Zhang { 6685582bec1SHong Zhang PetscErrorCode ierr; 6695582bec1SHong Zhang Mat_MLShell *shell; 6705582bec1SHong Zhang 6715582bec1SHong Zhang PetscFunctionBegin; 6725582bec1SHong Zhang ierr = MatShellGetContext(A,(void *)&shell);CHKERRQ(ierr); 6735582bec1SHong Zhang ierr = VecDestroy(shell->y);CHKERRQ(ierr); 6745582bec1SHong Zhang ierr = PetscFree(shell);CHKERRQ(ierr); 6755582bec1SHong Zhang ierr = MatDestroy_Shell(A);CHKERRQ(ierr); 6765582bec1SHong Zhang PetscFunctionReturn(0); 6775582bec1SHong Zhang } 6785582bec1SHong Zhang 6795582bec1SHong Zhang #undef __FUNCT__ 6805582bec1SHong Zhang #define __FUNCT__ "MatConvert_ML_SeqAIJ" 6815582bec1SHong Zhang PetscErrorCode MatConvert_ML_SeqAIJ(FineGridCtx *ml,Mat *newmat) 6825582bec1SHong Zhang { 6835582bec1SHong Zhang PetscErrorCode ierr; 6845582bec1SHong Zhang PetscInt i; 6855582bec1SHong Zhang ML_Operator *mat=ml->mlmat; 6865582bec1SHong Zhang PetscInt m=mat->outvec_leng,n= mat->invec_leng,nnz[m]; 6875582bec1SHong Zhang 6885582bec1SHong Zhang PetscFunctionBegin; 6895582bec1SHong Zhang if ( mat->getrow == NULL) SETERRQ(PETSC_ERR_ARG_NULL,"mat->getrow = NULL"); 690*9cb0ec93SHong Zhang #ifdef TEST 691*9cb0ec93SHong Zhang /* ---- new --------- */ 692*9cb0ec93SHong Zhang PetscInt *ai,*aj; 693*9cb0ec93SHong Zhang PetscScalar *aa; 694*9cb0ec93SHong Zhang ML_Matrix_DCSR *matdata = (ML_Matrix_DCSR*)mat->data; 6955582bec1SHong Zhang 696*9cb0ec93SHong Zhang ai = matdata->mat_ia; 697*9cb0ec93SHong Zhang aj = matdata->mat_ja; 698*9cb0ec93SHong Zhang aa = matdata->mat_a; 699*9cb0ec93SHong Zhang 700*9cb0ec93SHong Zhang /* -------- endof new ---------*/ 701*9cb0ec93SHong Zhang #endif /* TEST */ 7025582bec1SHong Zhang ierr = MatCreate(PETSC_COMM_SELF,m,n,PETSC_DECIDE,PETSC_DECIDE,newmat);CHKERRQ(ierr); 7035582bec1SHong Zhang ierr = MatSetType(*newmat,MATSEQAIJ);CHKERRQ(ierr); 7045582bec1SHong Zhang for (i=0; i<m; i++){ 7055582bec1SHong Zhang ML_get_matrix_row(mat,1,&i,&ml->rlen_max,&ml->cols,&ml->vals,&nnz[i],0); 7065582bec1SHong Zhang } 7075582bec1SHong Zhang ierr = MatSeqAIJSetPreallocation(*newmat,0,nnz);CHKERRQ(ierr); 7085582bec1SHong Zhang 7095582bec1SHong Zhang for (i=0; i<m; i++){ 7105582bec1SHong Zhang ML_get_matrix_row(mat,1,&i,&ml->rlen_max,&ml->cols,&ml->vals,&nnz[i],0); 7115582bec1SHong Zhang ierr = MatSetValues(*newmat,1,&i,nnz[i],ml->cols,ml->vals,INSERT_VALUES);CHKERRQ(ierr); 7125582bec1SHong Zhang } 7135582bec1SHong Zhang ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 7145582bec1SHong Zhang ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 7155582bec1SHong Zhang PetscFunctionReturn(0); 7165582bec1SHong Zhang } 7175582bec1SHong Zhang 7185582bec1SHong Zhang #undef __FUNCT__ 7195582bec1SHong Zhang #define __FUNCT__ "MatConvert_ML_SHELL" 7205582bec1SHong Zhang PetscErrorCode MatConvert_ML_SHELL(ML_Operator *mlmat,Mat *newmat) 7215582bec1SHong Zhang { 7225582bec1SHong Zhang PetscErrorCode ierr; 7235582bec1SHong Zhang PetscInt m,n; 7245582bec1SHong Zhang ML_Comm *MLcomm; 7255582bec1SHong Zhang Mat_MLShell *shellctx; 7265582bec1SHong Zhang 7275582bec1SHong Zhang PetscFunctionBegin; 7285582bec1SHong Zhang m = mlmat->outvec_leng; 7295582bec1SHong Zhang n = mlmat->invec_leng; 7305582bec1SHong Zhang if (!m || !n){ 7315582bec1SHong Zhang newmat = PETSC_NULL; 7325582bec1SHong Zhang } else { 7335582bec1SHong Zhang MLcomm = mlmat->comm; 7345582bec1SHong Zhang ierr = PetscNew(Mat_MLShell,&shellctx);CHKERRQ(ierr); 7355582bec1SHong Zhang ierr = PetscMemzero(shellctx,sizeof(Mat_MLShell));CHKERRQ(ierr); 7365582bec1SHong Zhang ierr = MatCreateShell(MLcomm->USR_comm,m,n,PETSC_DETERMINE,PETSC_DETERMINE,shellctx,newmat);CHKERRQ(ierr); 7375582bec1SHong Zhang ierr = MatShellSetOperation(*newmat,MATOP_MULT,(void *)MatMult_ML);CHKERRQ(ierr); 7385582bec1SHong Zhang ierr = MatShellSetOperation(*newmat,MATOP_MULT_ADD,(void *)MatMultAdd_ML);CHKERRQ(ierr); 7395582bec1SHong Zhang shellctx->A = *newmat; 7405582bec1SHong Zhang shellctx->mlmat = mlmat; 7415582bec1SHong Zhang ierr = VecCreate(PETSC_COMM_WORLD,&shellctx->y);CHKERRQ(ierr); 7425582bec1SHong Zhang ierr = VecSetSizes(shellctx->y,m,PETSC_DECIDE);CHKERRQ(ierr); 7435582bec1SHong Zhang ierr = VecSetFromOptions(shellctx->y);CHKERRQ(ierr); 7445582bec1SHong Zhang (*newmat)->ops->destroy = MatDestroy_ML; 7455582bec1SHong Zhang } 7465582bec1SHong Zhang PetscFunctionReturn(0); 7475582bec1SHong Zhang } 748