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