1*ab718edeSHong Zhang 25582bec1SHong Zhang /* 35582bec1SHong Zhang Provides an interface to the ML 3.0 smoothed Aggregation 45582bec1SHong Zhang */ 55582bec1SHong Zhang #include "src/ksp/pc/pcimpl.h" /*I "petscpc.h" I*/ 65582bec1SHong Zhang #include "src/ksp/pc/impls/mg/mgimpl.h" /*I "petscmg.h" I*/ 75582bec1SHong Zhang #include "src/mat/impls/aij/seq/aij.h" 85582bec1SHong Zhang #include "src/mat/impls/aij/mpi/mpiaij.h" 95582bec1SHong Zhang EXTERN_C_BEGIN 105582bec1SHong Zhang #include <math.h> 115582bec1SHong Zhang #include "ml_include.h" 125582bec1SHong Zhang EXTERN_C_END 135582bec1SHong Zhang 145582bec1SHong Zhang /* The context (data structure) at each grid level */ 155582bec1SHong Zhang typedef struct { 165582bec1SHong Zhang Vec x,b,r; /* global vectors */ 175582bec1SHong Zhang Mat A,P,R; 185582bec1SHong Zhang KSP ksp; 195582bec1SHong Zhang } GridCtx; 205582bec1SHong Zhang 215582bec1SHong Zhang /* The context used to input PETSc matrix into ML at fine grid */ 225582bec1SHong Zhang typedef struct { 235582bec1SHong Zhang Mat A,Aloc; 2424a42b14SHong Zhang Vec x,y; 255582bec1SHong Zhang ML_Operator *mlmat; 265582bec1SHong Zhang PetscScalar *pwork; /* tmp array used by PetscML_comm() */ 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 PetscErrorCode (*PCSetUp)(PC); 475582bec1SHong Zhang PetscErrorCode (*PCDestroy)(PC); 485582bec1SHong Zhang } PC_ML; 495582bec1SHong Zhang extern int PetscML_getrow(void *ML_data,int N_requested_rows,int requested_rows[], 505582bec1SHong Zhang int allocated_space,int columns[],double values[],int row_lengths[]); 515582bec1SHong Zhang extern int PetscML_matvec(void *ML_data, int in_length, double p[], int out_length,double ap[]); 525582bec1SHong Zhang extern int PetscML_comm(double x[], void *ML_data); 535582bec1SHong Zhang extern PetscErrorCode MatMult_ML(Mat,Vec,Vec); 545582bec1SHong Zhang extern PetscErrorCode MatMultAdd_ML(Mat,Vec,Vec,Vec); 555582bec1SHong Zhang extern PetscErrorCode MatConvert_MPIAIJ_ML(Mat,const MatType,Mat*); 565582bec1SHong Zhang extern PetscErrorCode MatDestroy_ML(Mat); 57e14861a4SHong Zhang extern PetscErrorCode MatConvert_ML_SeqAIJ(ML_Operator*,Mat*); 58*ab718edeSHong Zhang extern PetscErrorCode MatConvert_ML_MPIAIJ(ML_Operator*,Mat*); 595582bec1SHong Zhang extern PetscErrorCode MatConvert_ML_SHELL(ML_Operator*,Mat*); 605582bec1SHong Zhang 615582bec1SHong Zhang /* -------------------------------------------------------------------------- */ 625582bec1SHong Zhang /* 635582bec1SHong Zhang PCSetUp_ML - Prepares for the use of the ML preconditioner 645582bec1SHong Zhang by setting data structures and options. 655582bec1SHong Zhang 665582bec1SHong Zhang Input Parameter: 675582bec1SHong Zhang . pc - the preconditioner context 685582bec1SHong Zhang 695582bec1SHong Zhang Application Interface Routine: PCSetUp() 705582bec1SHong Zhang 715582bec1SHong Zhang Notes: 725582bec1SHong Zhang The interface routine PCSetUp() is not usually called directly by 735582bec1SHong Zhang the user, but instead is called by PCApply() if necessary. 745582bec1SHong Zhang */ 755582bec1SHong Zhang 765582bec1SHong Zhang #undef __FUNCT__ 775582bec1SHong Zhang #define __FUNCT__ "PCSetUp_ML" 785582bec1SHong Zhang static PetscErrorCode PCSetUp_ML(PC pc) 795582bec1SHong Zhang { 805582bec1SHong Zhang PetscErrorCode ierr; 815582bec1SHong Zhang PetscMPIInt size,rank; 825582bec1SHong Zhang FineGridCtx *PetscMLdata; 835582bec1SHong Zhang ML *ml_object; 845582bec1SHong Zhang ML_Aggregate *agg_object; 855582bec1SHong Zhang ML_Operator *mlmat; 865582bec1SHong Zhang PetscInt nlocal_allcols,Nlevels,mllevel,level,m,fine_level; 875582bec1SHong Zhang Mat A,Aloc; 885582bec1SHong Zhang GridCtx *gridctx; 895582bec1SHong Zhang PC_ML *pc_ml=PETSC_NULL; 905582bec1SHong Zhang PetscObjectContainer container; 915582bec1SHong Zhang 925582bec1SHong Zhang PetscFunctionBegin; 935582bec1SHong Zhang ierr = PetscObjectQuery((PetscObject)pc,"PC_ML",(PetscObject *)&container);CHKERRQ(ierr); 945582bec1SHong Zhang if (container) { 955582bec1SHong Zhang ierr = PetscObjectContainerGetPointer(container,(void **)&pc_ml);CHKERRQ(ierr); 965582bec1SHong Zhang } else { 975582bec1SHong Zhang SETERRQ(PETSC_ERR_ARG_NULL,"Container does not exit"); 985582bec1SHong Zhang } 995582bec1SHong Zhang 1005582bec1SHong Zhang /* setup special features of PCML */ 1015582bec1SHong Zhang /*--------------------------------*/ 1025582bec1SHong Zhang /* covert A to Aloc to be used by ML at fine grid */ 1035582bec1SHong Zhang A = pc->pmat; 1045582bec1SHong Zhang ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr); 105e14861a4SHong Zhang ierr = MPI_Comm_rank(A->comm,&rank);CHKERRQ(ierr); /* rm! */ 1065582bec1SHong Zhang pc_ml->size = size; 1075582bec1SHong Zhang if (size > 1){ 1085582bec1SHong Zhang Aloc = PETSC_NULL; 1095582bec1SHong Zhang ierr = MatConvert_MPIAIJ_ML(A,0,&Aloc);CHKERRQ(ierr); 1105582bec1SHong Zhang } else { 1115582bec1SHong Zhang Aloc = A; 1125582bec1SHong Zhang } 1135582bec1SHong Zhang 1145582bec1SHong Zhang /* create and initialize struct 'PetscMLdata' */ 1155582bec1SHong Zhang ierr = PetscNew(FineGridCtx,&PetscMLdata);CHKERRQ(ierr); 1165582bec1SHong Zhang ierr = PetscMemzero(PetscMLdata,sizeof(FineGridCtx));CHKERRQ(ierr); 1175582bec1SHong Zhang PetscMLdata->A = A; 1185582bec1SHong Zhang PetscMLdata->Aloc = Aloc; 1195582bec1SHong Zhang ierr = PetscMalloc((Aloc->n+1)*sizeof(PetscScalar),&PetscMLdata->pwork);CHKERRQ(ierr); 1205582bec1SHong Zhang pc_ml->PetscMLdata = PetscMLdata; 1215582bec1SHong Zhang 12224a42b14SHong Zhang ierr = VecCreate(PETSC_COMM_SELF,&PetscMLdata->x);CHKERRQ(ierr); 12324a42b14SHong Zhang if (size == 1){ 12424a42b14SHong Zhang ierr = VecSetSizes(PetscMLdata->x,A->n,PETSC_DECIDE);CHKERRQ(ierr); 12524a42b14SHong Zhang } else { 12624a42b14SHong Zhang ierr = VecSetSizes(PetscMLdata->x,Aloc->n,PETSC_DECIDE);CHKERRQ(ierr); 12724a42b14SHong Zhang } 12824a42b14SHong Zhang ierr = VecSetType(PetscMLdata->x,VECSEQ);CHKERRQ(ierr); 12924a42b14SHong Zhang 13024a42b14SHong Zhang ierr = VecCreate(PETSC_COMM_SELF,&PetscMLdata->y);CHKERRQ(ierr); 13124a42b14SHong Zhang ierr = VecSetSizes(PetscMLdata->y,A->m,PETSC_DECIDE);CHKERRQ(ierr); 13224a42b14SHong Zhang ierr = VecSetType(PetscMLdata->y,VECSEQ);CHKERRQ(ierr); 13324a42b14SHong Zhang 1345582bec1SHong Zhang /* create ML discretization matrix at fine grid */ 1355582bec1SHong Zhang ierr = MatGetSize(Aloc,&m,&nlocal_allcols);CHKERRQ(ierr); 1365582bec1SHong Zhang ML_Create(&ml_object,pc_ml->MaxNlevels); 1375582bec1SHong Zhang ML_Init_Amatrix(ml_object,0,m,m,PetscMLdata); 1385582bec1SHong Zhang ML_Set_Amatrix_Getrow(ml_object,0,PetscML_getrow,PetscML_comm,nlocal_allcols); 1395582bec1SHong Zhang ML_Set_Amatrix_Matvec(ml_object,0,PetscML_matvec); 1405582bec1SHong Zhang 1415582bec1SHong Zhang /* aggregation */ 1425582bec1SHong Zhang ML_Aggregate_Create(&agg_object); 1435582bec1SHong Zhang ML_Aggregate_Set_MaxCoarseSize(agg_object,pc_ml->MaxCoarseSize); 1445582bec1SHong Zhang /* set options */ 1455582bec1SHong Zhang switch (pc_ml->CoarsenScheme) { 1465582bec1SHong Zhang case 1: 1475582bec1SHong Zhang ML_Aggregate_Set_CoarsenScheme_Coupled(agg_object);break; 1485582bec1SHong Zhang case 2: 1495582bec1SHong Zhang ML_Aggregate_Set_CoarsenScheme_MIS(agg_object);break; 1505582bec1SHong Zhang case 3: 1515582bec1SHong Zhang ML_Aggregate_Set_CoarsenScheme_METIS(agg_object);break; 1525582bec1SHong Zhang } 1535582bec1SHong Zhang ML_Aggregate_Set_Threshold(agg_object,pc_ml->Threshold); 1545582bec1SHong Zhang ML_Aggregate_Set_DampingFactor(agg_object,pc_ml->DampingFactor); 1555582bec1SHong Zhang if (pc_ml->SpectralNormScheme_Anorm){ 1565582bec1SHong Zhang ML_Aggregate_Set_SpectralNormScheme_Anorm(agg_object); 1575582bec1SHong Zhang } 1585582bec1SHong Zhang 1595582bec1SHong Zhang Nlevels = ML_Gen_MGHierarchy_UsingAggregation(ml_object,0,ML_INCREASING,agg_object); 1605582bec1SHong Zhang if (Nlevels<=0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Nlevels %d must > 0",Nlevels); 1615582bec1SHong Zhang ierr = MGSetLevels(pc,Nlevels,PETSC_NULL);CHKERRQ(ierr); 1625582bec1SHong Zhang pc_ml->ml_object = ml_object; 1635582bec1SHong Zhang pc_ml->agg_object = agg_object; 1645582bec1SHong Zhang 1655582bec1SHong Zhang ierr = PetscMalloc(Nlevels*sizeof(GridCtx),&gridctx);CHKERRQ(ierr); 1665582bec1SHong Zhang fine_level = Nlevels - 1; 1675582bec1SHong Zhang pc_ml->gridctx = gridctx; 1685582bec1SHong Zhang pc_ml->fine_level = fine_level; 1695582bec1SHong Zhang 1705582bec1SHong Zhang /* wrap ML matrices by PETSc shell matrices at coarsened grids. 1715582bec1SHong Zhang Level 0 is the finest grid for ML, but coarsest for PETSc! */ 172e14861a4SHong Zhang gridctx[fine_level].A = A; 173e14861a4SHong Zhang level = fine_level - 1; 174*ab718edeSHong Zhang if (size == 1){ /* convert ML P, R and A into seqaij format */ 1755582bec1SHong Zhang for (mllevel=1; mllevel<Nlevels; mllevel++){ 176e14861a4SHong Zhang mlmat = &(ml_object->Pmat[mllevel]); 177e14861a4SHong Zhang ierr = MatConvert_ML_SeqAIJ(mlmat,&gridctx[level].P);CHKERRQ(ierr); 178e14861a4SHong Zhang mlmat = &(ml_object->Amat[mllevel]); 179e14861a4SHong Zhang ierr = MatConvert_ML_SeqAIJ(mlmat,&gridctx[level].A);CHKERRQ(ierr); 180e14861a4SHong Zhang mlmat = &(ml_object->Rmat[mllevel-1]); 181e14861a4SHong Zhang ierr = MatConvert_ML_SeqAIJ(mlmat,&gridctx[level].R);CHKERRQ(ierr); 1825582bec1SHong Zhang level--; 1835582bec1SHong Zhang } 184*ab718edeSHong Zhang } else { /* convert ML P and R into shell format, ML A into mpiaij format */ 1855582bec1SHong Zhang for (mllevel=1; mllevel<Nlevels; mllevel++){ 1865582bec1SHong Zhang mlmat = &(ml_object->Pmat[mllevel]); 1875582bec1SHong Zhang ierr = MatConvert_ML_SHELL(mlmat,&gridctx[level].P);CHKERRQ(ierr); 188*ab718edeSHong Zhang mlmat = &(ml_object->Rmat[mllevel-1]); 189*ab718edeSHong Zhang ierr = MatConvert_ML_SHELL(mlmat,&gridctx[level].R);CHKERRQ(ierr); 1905582bec1SHong Zhang mlmat = &(ml_object->Amat[mllevel]); 191*ab718edeSHong Zhang ierr = MatConvert_ML_MPIAIJ(mlmat,&gridctx[level].A);CHKERRQ(ierr); 192*ab718edeSHong Zhang #ifdef TMP 193e14861a4SHong Zhang /* 194e14861a4SHong Zhang if (mllevel==1) { 195e14861a4SHong Zhang ML_Operator_Print_UsingGlobalOrdering(mlmat,"Amat1",PETSC_NULL,PETSC_NULL); 196e14861a4SHong Zhang } 197e14861a4SHong Zhang */ 1985582bec1SHong Zhang ierr = MatConvert_ML_SHELL(mlmat,&gridctx[level].A);CHKERRQ(ierr); 199e14861a4SHong Zhang if (mllevel>0){ 200e14861a4SHong Zhang PetscMLdata->mlmat = &(ml_object->Amat[mllevel]); 201e14861a4SHong Zhang Mat A_tmp; 202e14861a4SHong Zhang ierr = MatConvert_ML_MPIAIJ(PetscMLdata,&A_tmp);CHKERRQ(ierr); 203e14861a4SHong Zhang /* ierr = MatView(A_tmp,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); */ 204e14861a4SHong Zhang 205e14861a4SHong Zhang Vec x,yy1,yy2; 206e14861a4SHong Zhang PetscInt am,an; 207e14861a4SHong Zhang ierr = MatGetLocalSize(A_tmp,&am,&an);CHKERRQ(ierr); 208e14861a4SHong Zhang ierr = VecCreate(PETSC_COMM_WORLD,&x);CHKERRQ(ierr); 209e14861a4SHong Zhang ierr = VecSetSizes(x,an,PETSC_DECIDE);CHKERRQ(ierr); 210e14861a4SHong Zhang ierr = VecSetFromOptions(x);CHKERRQ(ierr); 211e14861a4SHong Zhang ierr = VecCreate(PETSC_COMM_WORLD,&yy1);CHKERRQ(ierr); 212e14861a4SHong Zhang ierr = VecSetSizes(yy1,am,PETSC_DECIDE);CHKERRQ(ierr); 213e14861a4SHong Zhang ierr = VecSetFromOptions(yy1);CHKERRQ(ierr); 214e14861a4SHong Zhang ierr = VecDuplicate(yy1,&yy2);CHKERRQ(ierr); 215e14861a4SHong Zhang 216e14861a4SHong Zhang PetscRandom rd; 217e14861a4SHong Zhang PetscReal rnorm; 218e14861a4SHong Zhang PetscScalar mone = -1.0; 219e14861a4SHong Zhang ierr = PetscRandomCreate(PETSC_COMM_WORLD,RANDOM_DEFAULT,&rd);CHKERRQ(ierr); 220e14861a4SHong Zhang ierr = VecSetRandom(rd,x);CHKERRQ(ierr); 221e14861a4SHong Zhang ierr = MatMult(gridctx[level].A,x,yy1);CHKERRQ(ierr); 222e14861a4SHong Zhang ierr = MatMult(A_tmp,x,yy2);CHKERRQ(ierr); 223e14861a4SHong Zhang ierr = VecAXPY(&mone,yy1,yy2);CHKERRQ(ierr); 224e14861a4SHong Zhang ierr = VecNorm(yy2,NORM_2,&rnorm);CHKERRQ(ierr); 225e14861a4SHong Zhang printf(" [%d] mllevel: %d rnorm %g\n",rank,mllevel,rnorm); 226e14861a4SHong Zhang 227e14861a4SHong Zhang ierr = MatDestroy(A_tmp);CHKERRQ(ierr); 228e14861a4SHong Zhang ierr = VecDestroy(x);CHKERRQ(ierr); 229e14861a4SHong Zhang ierr = VecDestroy(yy1);CHKERRQ(ierr); 230e14861a4SHong Zhang ierr = VecDestroy(yy2);CHKERRQ(ierr); 231e14861a4SHong Zhang ierr = PetscRandomDestroy(rd);CHKERRQ(ierr); 232e14861a4SHong Zhang } 233e14861a4SHong Zhang #endif 2345582bec1SHong Zhang level--; 2355582bec1SHong Zhang } 2365582bec1SHong Zhang } 2375582bec1SHong Zhang 2385582bec1SHong Zhang /* create coarse level and the interpolation between the levels */ 2395582bec1SHong Zhang level = fine_level; 2405582bec1SHong Zhang while ( level >= 0 ){ 2415582bec1SHong Zhang if (level != fine_level){ 2425582bec1SHong Zhang ierr = VecCreate(gridctx[level].A->comm,&gridctx[level].x);CHKERRQ(ierr); 2435582bec1SHong Zhang ierr = VecSetSizes(gridctx[level].x,gridctx[level].A->n,PETSC_DECIDE);CHKERRQ(ierr); 2445582bec1SHong Zhang ierr = VecSetType(gridctx[level].x,VECMPI);CHKERRQ(ierr); 2455582bec1SHong Zhang ierr = MGSetX(pc,level,gridctx[level].x);CHKERRQ(ierr); 2465582bec1SHong Zhang 2475582bec1SHong Zhang ierr = VecCreate(gridctx[level].A->comm,&gridctx[level].b);CHKERRQ(ierr); 2485582bec1SHong Zhang ierr = VecSetSizes(gridctx[level].b,gridctx[level].A->m,PETSC_DECIDE);CHKERRQ(ierr); 2495582bec1SHong Zhang ierr = VecSetType(gridctx[level].b,VECMPI);CHKERRQ(ierr); 2505582bec1SHong Zhang ierr = MGSetRhs(pc,level,gridctx[level].b);CHKERRQ(ierr); 2515582bec1SHong Zhang } 2525582bec1SHong Zhang ierr = VecCreate(gridctx[level].A->comm,&gridctx[level].r);CHKERRQ(ierr); 2535582bec1SHong Zhang ierr = VecSetSizes(gridctx[level].r,gridctx[level].A->m,PETSC_DECIDE);CHKERRQ(ierr); 2545582bec1SHong Zhang ierr = VecSetType(gridctx[level].r,VECMPI);CHKERRQ(ierr); 2555582bec1SHong Zhang ierr = MGSetR(pc,level,gridctx[level].r);CHKERRQ(ierr); 2565582bec1SHong Zhang 2575582bec1SHong Zhang if (level == 0){ 2585582bec1SHong Zhang ierr = MGGetCoarseSolve(pc,&gridctx[level].ksp);CHKERRQ(ierr); 2595582bec1SHong Zhang } else { 2605582bec1SHong Zhang ierr = MGGetSmoother(pc,level,&gridctx[level].ksp);CHKERRQ(ierr); 2615582bec1SHong Zhang ierr = MGSetResidual(pc,level,MGDefaultResidual,gridctx[level].A);CHKERRQ(ierr); 2625582bec1SHong Zhang if (level == fine_level){ 2635582bec1SHong Zhang ierr = KSPSetOptionsPrefix(gridctx[level].ksp,"mg_fine_");CHKERRQ(ierr); 2645582bec1SHong Zhang ierr = MGSetR(pc,level,gridctx[level].r);CHKERRQ(ierr); 2655582bec1SHong Zhang } 2665582bec1SHong Zhang } 2675582bec1SHong Zhang ierr = KSPSetOperators(gridctx[level].ksp,gridctx[level].A,gridctx[level].A,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 2685582bec1SHong Zhang 2695582bec1SHong Zhang if (level < fine_level){ 2705582bec1SHong Zhang ierr = MGSetInterpolate(pc,level+1,gridctx[level].P);CHKERRQ(ierr); 2715582bec1SHong Zhang ierr = MGSetRestriction(pc,level+1,gridctx[level].R);CHKERRQ(ierr); 2725582bec1SHong Zhang } 2735582bec1SHong Zhang level--; 2745582bec1SHong Zhang } 2755582bec1SHong Zhang 2765582bec1SHong Zhang /* now call PCSetUp_MG() */ 2775582bec1SHong Zhang /*--------------------------------*/ 2785582bec1SHong Zhang ierr = (*pc_ml->PCSetUp)(pc);CHKERRQ(ierr); 2795582bec1SHong Zhang PetscFunctionReturn(0); 2805582bec1SHong Zhang } 2815582bec1SHong Zhang 2825582bec1SHong Zhang #undef __FUNCT__ 2835582bec1SHong Zhang #define __FUNCT__ "PetscObjectContainerDestroy_PC_ML" 2845582bec1SHong Zhang PetscErrorCode PetscObjectContainerDestroy_PC_ML(void *ptr) 2855582bec1SHong Zhang { 2865582bec1SHong Zhang PetscErrorCode ierr; 2875582bec1SHong Zhang PC_ML *pc_ml = (PC_ML*)ptr; 2885582bec1SHong Zhang PetscInt level; 2895582bec1SHong Zhang 2905582bec1SHong Zhang PetscFunctionBegin; 2915582bec1SHong Zhang if (pc_ml->size > 1){ierr = MatDestroy(pc_ml->PetscMLdata->Aloc);CHKERRQ(ierr);} 2925582bec1SHong Zhang ML_Aggregate_Destroy(&pc_ml->agg_object); 2935582bec1SHong Zhang ML_Destroy(&pc_ml->ml_object); 2945582bec1SHong Zhang 2955582bec1SHong Zhang ierr = PetscFree(pc_ml->PetscMLdata->pwork);CHKERRQ(ierr); 29624a42b14SHong Zhang ierr = VecDestroy(pc_ml->PetscMLdata->x);CHKERRQ(ierr); 29724a42b14SHong Zhang ierr = VecDestroy(pc_ml->PetscMLdata->y);CHKERRQ(ierr); 2985582bec1SHong Zhang ierr = PetscFree(pc_ml->PetscMLdata);CHKERRQ(ierr); 2995582bec1SHong Zhang 3005582bec1SHong Zhang level = pc_ml->fine_level; 3015582bec1SHong Zhang while ( level>= 0 ){ 3025582bec1SHong Zhang if (level != pc_ml->fine_level){ 3035582bec1SHong Zhang ierr = MatDestroy(pc_ml->gridctx[level].A);CHKERRQ(ierr); 3045582bec1SHong Zhang ierr = MatDestroy(pc_ml->gridctx[level].P);CHKERRQ(ierr); 3055582bec1SHong Zhang ierr = MatDestroy(pc_ml->gridctx[level].R);CHKERRQ(ierr); 3065582bec1SHong Zhang ierr = VecDestroy(pc_ml->gridctx[level].x);CHKERRQ(ierr); 3075582bec1SHong Zhang ierr = VecDestroy(pc_ml->gridctx[level].b);CHKERRQ(ierr); 3085582bec1SHong Zhang } 3095582bec1SHong Zhang ierr = VecDestroy(pc_ml->gridctx[level].r);CHKERRQ(ierr); 3105582bec1SHong Zhang level--; 3115582bec1SHong Zhang } 3125582bec1SHong Zhang ierr = PetscFree(pc_ml->gridctx);CHKERRQ(ierr); 3135582bec1SHong Zhang ierr = PetscFree(pc_ml);CHKERRQ(ierr); 3145582bec1SHong Zhang PetscFunctionReturn(0); 3155582bec1SHong Zhang } 3165582bec1SHong Zhang /* -------------------------------------------------------------------------- */ 3175582bec1SHong Zhang /* 3185582bec1SHong Zhang PCDestroy_ML - Destroys the private context for the ML preconditioner 3195582bec1SHong Zhang that was created with PCCreate_ML(). 3205582bec1SHong Zhang 3215582bec1SHong Zhang Input Parameter: 3225582bec1SHong Zhang . pc - the preconditioner context 3235582bec1SHong Zhang 3245582bec1SHong Zhang Application Interface Routine: PCDestroy() 3255582bec1SHong Zhang */ 3265582bec1SHong Zhang #undef __FUNCT__ 3275582bec1SHong Zhang #define __FUNCT__ "PCDestroy_ML" 3285582bec1SHong Zhang static PetscErrorCode PCDestroy_ML(PC pc) 3295582bec1SHong Zhang { 3305582bec1SHong Zhang PetscErrorCode ierr; 3315582bec1SHong Zhang PC_ML *pc_ml=PETSC_NULL; 3325582bec1SHong Zhang PetscObjectContainer container; 3335582bec1SHong Zhang 3345582bec1SHong Zhang PetscFunctionBegin; 3355582bec1SHong Zhang ierr = PetscObjectQuery((PetscObject)pc,"PC_ML",(PetscObject *)&container);CHKERRQ(ierr); 3365582bec1SHong Zhang if (container) { 3375582bec1SHong Zhang ierr = PetscObjectContainerGetPointer(container,(void **)&pc_ml);CHKERRQ(ierr); 3385582bec1SHong Zhang pc->ops->destroy = pc_ml->PCDestroy; 3395582bec1SHong Zhang } else { 3405582bec1SHong Zhang SETERRQ(PETSC_ERR_ARG_NULL,"Container does not exit"); 3415582bec1SHong Zhang } 3429cb0ec93SHong Zhang /* detach pc and PC_ML and dereference container */ 3435582bec1SHong Zhang ierr = PetscObjectCompose((PetscObject)pc,"PC_ML",0);CHKERRQ(ierr); 3445582bec1SHong Zhang ierr = (*pc->ops->destroy)(pc);CHKERRQ(ierr); 3455582bec1SHong Zhang 3465582bec1SHong Zhang ierr = PetscObjectContainerDestroy(container);CHKERRQ(ierr); 3475582bec1SHong Zhang PetscFunctionReturn(0); 3485582bec1SHong Zhang } 3495582bec1SHong Zhang 3505582bec1SHong Zhang #undef __FUNCT__ 3515582bec1SHong Zhang #define __FUNCT__ "PCSetFromOptions_ML" 3525582bec1SHong Zhang static PetscErrorCode PCSetFromOptions_ML(PC pc) 3535582bec1SHong Zhang { 3545582bec1SHong Zhang PetscErrorCode ierr; 3555582bec1SHong Zhang PetscInt indx,m,PrintLevel,MaxNlevels,MaxCoarseSize; 3565582bec1SHong Zhang PetscReal Threshold,DampingFactor; 3575582bec1SHong Zhang PetscTruth flg; 3585582bec1SHong Zhang const char *type[] = {"additive","multiplicative","full","cascade","kascade"}; 3595582bec1SHong Zhang const char *scheme[] = {"Uncoupled","Coupled","MIS","METIS"}; 3605582bec1SHong Zhang PC_ML *pc_ml=PETSC_NULL; 3615582bec1SHong Zhang PetscObjectContainer container; 3625582bec1SHong Zhang 3635582bec1SHong Zhang PetscFunctionBegin; 3645582bec1SHong Zhang ierr = PetscObjectQuery((PetscObject)pc,"PC_ML",(PetscObject *)&container);CHKERRQ(ierr); 3655582bec1SHong Zhang if (container) { 3665582bec1SHong Zhang ierr = PetscObjectContainerGetPointer(container,(void **)&pc_ml);CHKERRQ(ierr); 3675582bec1SHong Zhang } else { 3685582bec1SHong Zhang SETERRQ(PETSC_ERR_ARG_NULL,"Container does not exit"); 3695582bec1SHong Zhang } 3705582bec1SHong Zhang ierr = PetscOptionsHead("MG options");CHKERRQ(ierr); 3715582bec1SHong Zhang /* inherited MG options */ 3725582bec1SHong Zhang ierr = PetscOptionsInt("-pc_mg_cycles","1 for V cycle, 2 for W-cycle","MGSetCycles",1,&m,&flg);CHKERRQ(ierr); 3735582bec1SHong Zhang if (flg) { 3745582bec1SHong Zhang ierr = MGSetCycles(pc,m);CHKERRQ(ierr); 3755582bec1SHong Zhang } 3765582bec1SHong Zhang ierr = PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","MGSetNumberSmoothUp",1,&m,&flg);CHKERRQ(ierr); 3775582bec1SHong Zhang if (flg) { 3785582bec1SHong Zhang ierr = MGSetNumberSmoothUp(pc,m);CHKERRQ(ierr); 3795582bec1SHong Zhang } 3805582bec1SHong Zhang ierr = PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","MGSetNumberSmoothDown",1,&m,&flg);CHKERRQ(ierr); 3815582bec1SHong Zhang if (flg) { 3825582bec1SHong Zhang ierr = MGSetNumberSmoothDown(pc,m);CHKERRQ(ierr); 3835582bec1SHong Zhang } 3845582bec1SHong Zhang ierr = PetscOptionsEList("-pc_mg_type","Multigrid type","MGSetType",type,5,type[1],&indx,&flg);CHKERRQ(ierr); 3855582bec1SHong Zhang if (flg) { 3865582bec1SHong Zhang MGType mg = (MGType) 0; 3875582bec1SHong Zhang switch (indx) { 3885582bec1SHong Zhang case 0: 3895582bec1SHong Zhang mg = MGADDITIVE; 3905582bec1SHong Zhang break; 3915582bec1SHong Zhang case 1: 3925582bec1SHong Zhang mg = MGMULTIPLICATIVE; 3935582bec1SHong Zhang break; 3945582bec1SHong Zhang case 2: 3955582bec1SHong Zhang mg = MGFULL; 3965582bec1SHong Zhang break; 3975582bec1SHong Zhang case 3: 3985582bec1SHong Zhang mg = MGKASKADE; 3995582bec1SHong Zhang break; 4005582bec1SHong Zhang case 4: 4015582bec1SHong Zhang mg = MGKASKADE; 4025582bec1SHong Zhang break; 4035582bec1SHong Zhang } 4045582bec1SHong Zhang ierr = MGSetType(pc,mg);CHKERRQ(ierr); 4055582bec1SHong Zhang } 4065582bec1SHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 4075582bec1SHong Zhang 4085582bec1SHong Zhang /* ML options */ 4095582bec1SHong Zhang ierr = PetscOptionsHead("ML options");CHKERRQ(ierr); 4105582bec1SHong Zhang /* set defaults */ 4115582bec1SHong Zhang PrintLevel = 0; 4125582bec1SHong Zhang MaxNlevels = 10; 4135582bec1SHong Zhang MaxCoarseSize = 1; 4145582bec1SHong Zhang indx = 0; 4155582bec1SHong Zhang Threshold = 0.0; 4165582bec1SHong Zhang DampingFactor = 4.0/3.0; 4175582bec1SHong Zhang 4185582bec1SHong Zhang ierr = PetscOptionsInt("-pc_ml_PrintLevel","Print level","ML_Set_PrintLevel",PrintLevel,&PrintLevel,PETSC_NULL);CHKERRQ(ierr); 4195582bec1SHong Zhang ML_Set_PrintLevel(PrintLevel); 4205582bec1SHong Zhang 4215582bec1SHong Zhang ierr = PetscOptionsInt("-pc_ml_maxNlevels","Maximum number of levels","None",MaxNlevels,&MaxNlevels,PETSC_NULL);CHKERRQ(ierr); 4225582bec1SHong Zhang pc_ml->MaxNlevels = MaxNlevels; 4235582bec1SHong Zhang 4245582bec1SHong Zhang ierr = PetscOptionsInt("-pc_ml_maxCoarseSize","Maximum coarsest mesh size","ML_Aggregate_Set_MaxCoarseSize",MaxCoarseSize,&MaxCoarseSize,PETSC_NULL);CHKERRQ(ierr); 4255582bec1SHong Zhang pc_ml->MaxCoarseSize = MaxCoarseSize; 4265582bec1SHong Zhang 4275582bec1SHong Zhang ierr = PetscOptionsEList("-pc_ml_CoarsenScheme","Aggregate Coarsen Scheme","ML_Aggregate_Set_CoarsenScheme_*",scheme,4,scheme[0],&indx,PETSC_NULL);CHKERRQ(ierr); 4285582bec1SHong Zhang pc_ml->CoarsenScheme = indx; 4295582bec1SHong Zhang 4305582bec1SHong Zhang ierr = PetscOptionsReal("-pc_ml_DampingFactor","P damping factor","ML_Aggregate_Set_DampingFactor",DampingFactor,&DampingFactor,PETSC_NULL);CHKERRQ(ierr); 4315582bec1SHong Zhang pc_ml->DampingFactor = DampingFactor; 4325582bec1SHong Zhang 4335582bec1SHong Zhang ierr = PetscOptionsReal("-pc_ml_Threshold","Smoother drop tol","ML_Aggregate_Set_Threshold",Threshold,&Threshold,PETSC_NULL);CHKERRQ(ierr); 4345582bec1SHong Zhang pc_ml->Threshold = Threshold; 4355582bec1SHong Zhang 4365582bec1SHong 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); 4375582bec1SHong Zhang 4385582bec1SHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 4395582bec1SHong Zhang PetscFunctionReturn(0); 4405582bec1SHong Zhang } 4415582bec1SHong Zhang 4425582bec1SHong Zhang /* -------------------------------------------------------------------------- */ 4435582bec1SHong Zhang /* 4445582bec1SHong Zhang PCCreate_ML - Creates a ML preconditioner context, PC_ML, 4455582bec1SHong Zhang and sets this as the private data within the generic preconditioning 4465582bec1SHong Zhang context, PC, that was created within PCCreate(). 4475582bec1SHong Zhang 4485582bec1SHong Zhang Input Parameter: 4495582bec1SHong Zhang . pc - the preconditioner context 4505582bec1SHong Zhang 4515582bec1SHong Zhang Application Interface Routine: PCCreate() 4525582bec1SHong Zhang */ 4535582bec1SHong Zhang 4545582bec1SHong Zhang /*MC 4555582bec1SHong Zhang PCML - Use geometric multigrid preconditioning. This preconditioner requires you provide 4565582bec1SHong Zhang fine grid discretization matrix. The coarser grid matrices and restriction/interpolation 4575582bec1SHong Zhang operators are computed by ML and wrapped as PETSc shell matrices. 4585582bec1SHong Zhang 4595582bec1SHong Zhang Options Database Key: (not done yet!) 4605582bec1SHong Zhang + -pc_mg_maxlevels <nlevels> - maximum number of levels including finest 4615582bec1SHong Zhang . -pc_mg_cycles 1 or 2 - for V or W-cycle 4625582bec1SHong Zhang . -pc_mg_smoothup <n> - number of smoothing steps after interpolation 4635582bec1SHong Zhang . -pc_mg_smoothdown <n> - number of smoothing steps before applying restriction operator 4645582bec1SHong Zhang . -pc_mg_type <additive,multiplicative,full,cascade> - multiplicative is the default 4655582bec1SHong Zhang . -pc_mg_monitor - print information on the multigrid convergence 4665582bec1SHong Zhang - -pc_mg_dump_matlab - dumps the matrices for each level and the restriction/interpolation matrices 4675582bec1SHong Zhang to the Socket viewer for reading from Matlab. 4685582bec1SHong Zhang 4695582bec1SHong Zhang Level: intermediate 4705582bec1SHong Zhang 4715582bec1SHong Zhang Concepts: multigrid 4725582bec1SHong Zhang 4735582bec1SHong Zhang .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, 4745582bec1SHong Zhang MGSetLevels(), MGGetLevels(), MGSetType(), MPSetCycles(), MGSetNumberSmoothDown(), 4755582bec1SHong Zhang MGSetNumberSmoothUp(), MGGetCoarseSolve(), MGSetResidual(), MGSetInterpolation(), 4765582bec1SHong Zhang MGSetRestriction(), MGGetSmoother(), MGGetSmootherUp(), MGGetSmootherDown(), 4775582bec1SHong Zhang MGSetCyclesOnLevel(), MGSetRhs(), MGSetX(), MGSetR() 4785582bec1SHong Zhang M*/ 4795582bec1SHong Zhang 4805582bec1SHong Zhang EXTERN_C_BEGIN 4815582bec1SHong Zhang #undef __FUNCT__ 4825582bec1SHong Zhang #define __FUNCT__ "PCCreate_ML" 4835582bec1SHong Zhang PetscErrorCode PCCreate_ML(PC pc) 4845582bec1SHong Zhang { 4855582bec1SHong Zhang PetscErrorCode ierr; 4865582bec1SHong Zhang PC_ML *pc_ml; 4875582bec1SHong Zhang PetscObjectContainer container; 4885582bec1SHong Zhang 4895582bec1SHong Zhang PetscFunctionBegin; 4905582bec1SHong Zhang /* initialize pc as PCMG */ 4915582bec1SHong Zhang ierr = PCSetType(pc,PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */ 4925582bec1SHong Zhang 4935582bec1SHong Zhang /* create a supporting struct and attach it to pc */ 4945582bec1SHong Zhang ierr = PetscNew(PC_ML,&pc_ml);CHKERRQ(ierr); 4959cb0ec93SHong Zhang ierr = PetscMemzero(pc_ml,sizeof(PC_ML));CHKERRQ(ierr); 4965582bec1SHong Zhang ierr = PetscObjectContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr); 4975582bec1SHong Zhang ierr = PetscObjectContainerSetPointer(container,pc_ml);CHKERRQ(ierr); 4989cb0ec93SHong Zhang ierr = PetscObjectContainerSetUserDestroy(container,PetscObjectContainerDestroy_PC_ML);CHKERRQ(ierr); 4995582bec1SHong Zhang ierr = PetscObjectCompose((PetscObject)pc,"PC_ML",(PetscObject)container);CHKERRQ(ierr); 5005582bec1SHong Zhang 5015582bec1SHong Zhang pc_ml->PCSetUp = pc->ops->setup; 5025582bec1SHong Zhang pc_ml->PCDestroy = pc->ops->destroy; 5035582bec1SHong Zhang 5045582bec1SHong Zhang /* overwrite the pointers of PCMG by the functions of PCML */ 5055582bec1SHong Zhang pc->ops->setfromoptions = PCSetFromOptions_ML; 5065582bec1SHong Zhang pc->ops->setup = PCSetUp_ML; 5075582bec1SHong Zhang pc->ops->destroy = PCDestroy_ML; 5085582bec1SHong Zhang PetscFunctionReturn(0); 5095582bec1SHong Zhang } 5105582bec1SHong Zhang EXTERN_C_END 5115582bec1SHong Zhang 5125582bec1SHong Zhang int PetscML_getrow(void *ML_data, int N_requested_rows, int requested_rows[], 5135582bec1SHong Zhang int allocated_space, int columns[], double values[], int row_lengths[]) 5145582bec1SHong Zhang { 5155582bec1SHong Zhang PetscErrorCode ierr; 5165582bec1SHong Zhang Mat Aloc; 5175582bec1SHong Zhang Mat_SeqAIJ *a; 5185582bec1SHong Zhang PetscInt m,i,j,k=0,row,*aj; 5195582bec1SHong Zhang PetscScalar *aa; 5205582bec1SHong Zhang FineGridCtx *ml=(FineGridCtx*)ML_data; 5215582bec1SHong Zhang 5225582bec1SHong Zhang Aloc = ml->Aloc; 5235582bec1SHong Zhang a = (Mat_SeqAIJ*)Aloc->data; 5245582bec1SHong Zhang ierr = MatGetSize(Aloc,&m,PETSC_NULL);CHKERRQ(ierr); 5255582bec1SHong Zhang 5265582bec1SHong Zhang for (i = 0; i<N_requested_rows; i++) { 5275582bec1SHong Zhang row = requested_rows[i]; 5285582bec1SHong Zhang row_lengths[i] = a->ilen[row]; 5295582bec1SHong Zhang if (allocated_space < k+row_lengths[i]) return(0); 5305582bec1SHong Zhang if ( (row >= 0) || (row <= (m-1)) ) { 5315582bec1SHong Zhang aj = a->j + a->i[row]; 5325582bec1SHong Zhang aa = a->a + a->i[row]; 5335582bec1SHong Zhang for (j=0; j<row_lengths[i]; j++){ 5345582bec1SHong Zhang columns[k] = aj[j]; 5355582bec1SHong Zhang values[k++] = aa[j]; 5365582bec1SHong Zhang } 5375582bec1SHong Zhang } 5385582bec1SHong Zhang } 5395582bec1SHong Zhang return(1); 5405582bec1SHong Zhang } 5415582bec1SHong Zhang 5425582bec1SHong Zhang int PetscML_matvec(void *ML_data,int in_length,double p[],int out_length,double ap[]) 5435582bec1SHong Zhang { 5445582bec1SHong Zhang PetscErrorCode ierr; 5455582bec1SHong Zhang FineGridCtx *ml=(FineGridCtx*)ML_data; 5465582bec1SHong Zhang Mat A=ml->A, Aloc=ml->Aloc; 5475582bec1SHong Zhang PetscMPIInt size; 5485582bec1SHong Zhang PetscScalar *pwork=ml->pwork; 5495582bec1SHong Zhang PetscInt i; 5505582bec1SHong Zhang 5515582bec1SHong Zhang ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr); 5525582bec1SHong Zhang if (size == 1){ 55324a42b14SHong Zhang ierr = VecPlaceArray(ml->x,p);CHKERRQ(ierr); 5545582bec1SHong Zhang } else { 5555582bec1SHong Zhang for (i=0; i<in_length; i++) pwork[i] = p[i]; 5565582bec1SHong Zhang PetscML_comm(pwork,ml); 55724a42b14SHong Zhang ierr = VecPlaceArray(ml->x,pwork);CHKERRQ(ierr); 5585582bec1SHong Zhang } 55924a42b14SHong Zhang ierr = VecPlaceArray(ml->y,ap);CHKERRQ(ierr); 56024a42b14SHong Zhang ierr = MatMult(Aloc,ml->x,ml->y);CHKERRQ(ierr); 5615582bec1SHong Zhang return 0; 5625582bec1SHong Zhang } 5635582bec1SHong Zhang 5645582bec1SHong Zhang int PetscML_comm(double p[],void *ML_data) 5655582bec1SHong Zhang { 5665582bec1SHong Zhang PetscErrorCode ierr; 5675582bec1SHong Zhang FineGridCtx *ml=(FineGridCtx*)ML_data; 5685582bec1SHong Zhang Mat A=ml->A; 5695582bec1SHong Zhang Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data; 5705582bec1SHong Zhang PetscMPIInt size; 5715582bec1SHong Zhang PetscInt i,in_length=A->m,out_length=ml->Aloc->n; 5725582bec1SHong Zhang PetscScalar *array; 5735582bec1SHong Zhang 5745582bec1SHong Zhang ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr); 5755582bec1SHong Zhang if (size == 1) return 0; 57624a42b14SHong Zhang 57724a42b14SHong Zhang ierr = VecPlaceArray(ml->y,p);CHKERRQ(ierr); 57824a42b14SHong Zhang ierr = VecScatterBegin(ml->y,a->lvec,INSERT_VALUES,SCATTER_FORWARD,a->Mvctx);CHKERRQ(ierr); 57924a42b14SHong Zhang ierr = VecScatterEnd(ml->y,a->lvec,INSERT_VALUES,SCATTER_FORWARD,a->Mvctx);CHKERRQ(ierr); 5805582bec1SHong Zhang ierr = VecGetArray(a->lvec,&array);CHKERRQ(ierr); 5815582bec1SHong Zhang for (i=in_length; i<out_length; i++){ 5825582bec1SHong Zhang p[i] = array[i-in_length]; 5835582bec1SHong Zhang } 5845582bec1SHong Zhang return 0; 5855582bec1SHong Zhang } 5865582bec1SHong Zhang #undef __FUNCT__ 5875582bec1SHong Zhang #define __FUNCT__ "MatMult_ML" 5885582bec1SHong Zhang PetscErrorCode MatMult_ML(Mat A,Vec x,Vec y) 5895582bec1SHong Zhang { 5905582bec1SHong Zhang PetscErrorCode ierr; 5915582bec1SHong Zhang Mat_MLShell *shell; 5925582bec1SHong Zhang PetscScalar *xarray,*yarray; 5935582bec1SHong Zhang PetscInt x_length,y_length; 5945582bec1SHong Zhang 5955582bec1SHong Zhang PetscFunctionBegin; 5965582bec1SHong Zhang ierr = MatShellGetContext(A,(void *)&shell);CHKERRQ(ierr); 5975582bec1SHong Zhang ierr = VecGetArray(x,&xarray);CHKERRQ(ierr); 5985582bec1SHong Zhang ierr = VecGetArray(y,&yarray);CHKERRQ(ierr); 5995582bec1SHong Zhang x_length = shell->mlmat->invec_leng; 6005582bec1SHong Zhang y_length = shell->mlmat->outvec_leng; 6015582bec1SHong Zhang 6025582bec1SHong Zhang ML_Operator_Apply(shell->mlmat,x_length,xarray,y_length,yarray); 6035582bec1SHong Zhang 6045582bec1SHong Zhang ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr); 6055582bec1SHong Zhang ierr = VecRestoreArray(y,&yarray);CHKERRQ(ierr); 6065582bec1SHong Zhang PetscFunctionReturn(0); 6075582bec1SHong Zhang } 6085582bec1SHong Zhang /* MatMultAdd_ML - Compute y = w + A*x */ 6095582bec1SHong Zhang #undef __FUNCT__ 6105582bec1SHong Zhang #define __FUNCT__ "MatMultAdd_ML" 6115582bec1SHong Zhang PetscErrorCode MatMultAdd_ML(Mat A,Vec x,Vec w,Vec y) 6125582bec1SHong Zhang { 6135582bec1SHong Zhang PetscErrorCode ierr; 6145582bec1SHong Zhang Mat_MLShell *shell; 6155582bec1SHong Zhang PetscScalar *xarray,*yarray; 6165582bec1SHong Zhang const PetscScalar one=1.0; 6175582bec1SHong Zhang PetscInt x_length,y_length; 6185582bec1SHong Zhang 6195582bec1SHong Zhang PetscFunctionBegin; 6205582bec1SHong Zhang ierr = MatShellGetContext(A,(void *)&shell);CHKERRQ(ierr); 6215582bec1SHong Zhang ierr = VecGetArray(x,&xarray);CHKERRQ(ierr); 6225582bec1SHong Zhang ierr = VecGetArray(y,&yarray);CHKERRQ(ierr); 6235582bec1SHong Zhang 6245582bec1SHong Zhang x_length = shell->mlmat->invec_leng; 6255582bec1SHong Zhang y_length = shell->mlmat->outvec_leng; 6265582bec1SHong Zhang 6275582bec1SHong Zhang ML_Operator_Apply(shell->mlmat,x_length,xarray,y_length,yarray); 6285582bec1SHong Zhang 6295582bec1SHong Zhang ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr); 6305582bec1SHong Zhang ierr = VecRestoreArray(y,&yarray);CHKERRQ(ierr); 6315582bec1SHong Zhang ierr = VecAXPY(&one,w,y);CHKERRQ(ierr); 6325582bec1SHong Zhang 6335582bec1SHong Zhang PetscFunctionReturn(0); 6345582bec1SHong Zhang } 6355582bec1SHong Zhang 6365582bec1SHong Zhang /* newtype is ignored because "ml" is not listed under Petsc MatType yet */ 6375582bec1SHong Zhang #undef __FUNCT__ 6385582bec1SHong Zhang #define __FUNCT__ "MatConvert_MPIAIJ_ML" 6395582bec1SHong Zhang PetscErrorCode MatConvert_MPIAIJ_ML(Mat A,const MatType newtype,Mat *Aloc) 6405582bec1SHong Zhang { 6415582bec1SHong Zhang PetscErrorCode ierr; 6425582bec1SHong Zhang Mat_MPIAIJ *mpimat=(Mat_MPIAIJ*)A->data; 6435582bec1SHong Zhang Mat_SeqAIJ *mat,*a=(Mat_SeqAIJ*)(mpimat->A)->data,*b=(Mat_SeqAIJ*)(mpimat->B)->data; 6445582bec1SHong Zhang PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j; 6455582bec1SHong Zhang PetscScalar *aa=a->a,*ba=b->a,*ca; 6465582bec1SHong Zhang PetscInt am=A->m,an=A->n,i,j,k; 6475582bec1SHong Zhang PetscInt *ci,*cj,ncols; 6485582bec1SHong Zhang MatReuse scall=MAT_INITIAL_MATRIX; 6495582bec1SHong Zhang 6505582bec1SHong Zhang PetscFunctionBegin; 6515582bec1SHong Zhang if (am != an) SETERRQ2(PETSC_ERR_ARG_WRONG,"A must have a square diagonal portion, am: %d != an: %d",am,an); 6525582bec1SHong Zhang 6535582bec1SHong Zhang if (*Aloc) scall = MAT_REUSE_MATRIX; 6545582bec1SHong Zhang if (scall == MAT_INITIAL_MATRIX){ 6555582bec1SHong Zhang ierr = PetscMalloc((1+am)*sizeof(PetscInt),&ci);CHKERRQ(ierr); 6565582bec1SHong Zhang ci[0] = 0; 6575582bec1SHong Zhang for (i=0; i<am; i++){ 6585582bec1SHong Zhang ci[i+1] = ci[i] + (ai[i+1] - ai[i]) + (bi[i+1] - bi[i]); 6595582bec1SHong Zhang } 6605582bec1SHong Zhang ierr = PetscMalloc((1+ci[am])*sizeof(PetscInt),&cj);CHKERRQ(ierr); 6615582bec1SHong Zhang ierr = PetscMalloc((1+ci[am])*sizeof(PetscScalar),&ca);CHKERRQ(ierr); 6625582bec1SHong Zhang 6635582bec1SHong Zhang k = 0; 6645582bec1SHong Zhang for (i=0; i<am; i++){ 6655582bec1SHong Zhang /* diagonal portion of A */ 6665582bec1SHong Zhang ncols = ai[i+1] - ai[i]; 6675582bec1SHong Zhang for (j=0; j<ncols; j++) { 6685582bec1SHong Zhang cj[k] = *aj++; 6695582bec1SHong Zhang ca[k++] = *aa++; 6705582bec1SHong Zhang } 6715582bec1SHong Zhang /* off-diagonal portion of A */ 6725582bec1SHong Zhang ncols = bi[i+1] - bi[i]; 6735582bec1SHong Zhang for (j=0; j<ncols; j++) { 6745582bec1SHong Zhang cj[k] = an + (*bj); bj++; 6755582bec1SHong Zhang ca[k++] = *ba++; 6765582bec1SHong Zhang } 6775582bec1SHong Zhang } 6785582bec1SHong Zhang if (k != ci[am]) SETERRQ2(PETSC_ERR_ARG_WRONG,"k: %d != ci[am]: %d",k,ci[am]); 6795582bec1SHong Zhang 6805582bec1SHong Zhang /* put together the new matrix */ 6815582bec1SHong Zhang an = mpimat->A->n+mpimat->B->n; 6825582bec1SHong Zhang ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,an,ci,cj,ca,Aloc);CHKERRQ(ierr); 6835582bec1SHong Zhang 6845582bec1SHong Zhang /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 6855582bec1SHong Zhang /* Since these are PETSc arrays, change flags to free them as necessary. */ 6865582bec1SHong Zhang mat = (Mat_SeqAIJ*)(*Aloc)->data; 6875582bec1SHong Zhang mat->freedata = PETSC_TRUE; 6885582bec1SHong Zhang mat->nonew = 0; 6895582bec1SHong Zhang } else if (scall == MAT_REUSE_MATRIX){ 6905582bec1SHong Zhang mat=(Mat_SeqAIJ*)(*Aloc)->data; 6915582bec1SHong Zhang ci = mat->i; cj = mat->j; ca = mat->a; 6925582bec1SHong Zhang for (i=0; i<am; i++) { 6935582bec1SHong Zhang /* diagonal portion of A */ 6945582bec1SHong Zhang ncols = ai[i+1] - ai[i]; 6955582bec1SHong Zhang for (j=0; j<ncols; j++) *ca++ = *aa++; 6965582bec1SHong Zhang /* off-diagonal portion of A */ 6975582bec1SHong Zhang ncols = bi[i+1] - bi[i]; 6985582bec1SHong Zhang for (j=0; j<ncols; j++) *ca++ = *ba++; 6995582bec1SHong Zhang } 7005582bec1SHong Zhang } else { 7015582bec1SHong Zhang SETERRQ1(PETSC_ERR_ARG_WRONG,"Invalid MatReuse %d",(int)scall); 7025582bec1SHong Zhang } 7035582bec1SHong Zhang PetscFunctionReturn(0); 7045582bec1SHong Zhang } 7055582bec1SHong Zhang extern PetscErrorCode MatDestroy_Shell(Mat); 7065582bec1SHong Zhang #undef __FUNCT__ 7075582bec1SHong Zhang #define __FUNCT__ "MatDestroy_ML" 7085582bec1SHong Zhang PetscErrorCode MatDestroy_ML(Mat A) 7095582bec1SHong Zhang { 7105582bec1SHong Zhang PetscErrorCode ierr; 7115582bec1SHong Zhang Mat_MLShell *shell; 7125582bec1SHong Zhang 7135582bec1SHong Zhang PetscFunctionBegin; 7145582bec1SHong Zhang ierr = MatShellGetContext(A,(void *)&shell);CHKERRQ(ierr); 7155582bec1SHong Zhang ierr = VecDestroy(shell->y);CHKERRQ(ierr); 7165582bec1SHong Zhang ierr = PetscFree(shell);CHKERRQ(ierr); 7175582bec1SHong Zhang ierr = MatDestroy_Shell(A);CHKERRQ(ierr); 7185582bec1SHong Zhang PetscFunctionReturn(0); 7195582bec1SHong Zhang } 7205582bec1SHong Zhang 721e14861a4SHong Zhang extern PetscErrorCode PetscSortIntWithScalarArray(PetscInt,PetscInt [],PetscScalar []); 7225582bec1SHong Zhang #undef __FUNCT__ 7235582bec1SHong Zhang #define __FUNCT__ "MatConvert_ML_SeqAIJ" 724e14861a4SHong Zhang PetscErrorCode MatConvert_ML_SeqAIJ(ML_Operator *mlmat,Mat *newmat) 7255582bec1SHong Zhang { 726e14861a4SHong Zhang struct ML_CSR_MSRdata *matdata = (struct ML_CSR_MSRdata *)mlmat->data; 7275582bec1SHong Zhang PetscErrorCode ierr; 728e14861a4SHong Zhang PetscInt m=mlmat->outvec_leng,n=mlmat->invec_leng,nnz[m],nz_max; 729e14861a4SHong Zhang PetscInt *ml_cols=matdata->columns,*aj,i,j,k; 730e14861a4SHong Zhang PetscScalar *ml_vals=matdata->values,*aa; 7315582bec1SHong Zhang 7325582bec1SHong Zhang PetscFunctionBegin; 733e14861a4SHong Zhang if ( mlmat->getrow == NULL) SETERRQ(PETSC_ERR_ARG_NULL,"mlmat->getrow = NULL"); 734*ab718edeSHong Zhang if (m != n){ /* ML Pmat and Rmat are in CSR format. Pass array pointers into SeqAIJ matrix */ 735e14861a4SHong Zhang ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,m,n,matdata->rowptr,ml_cols,ml_vals,newmat);CHKERRQ(ierr); 73624a42b14SHong Zhang PetscFunctionReturn(0); 73724a42b14SHong Zhang } 7385582bec1SHong Zhang 739e14861a4SHong Zhang /* ML Amat is in MSR format. Copy its data into SeqAIJ matrix */ 7405582bec1SHong Zhang ierr = MatCreate(PETSC_COMM_SELF,m,n,PETSC_DECIDE,PETSC_DECIDE,newmat);CHKERRQ(ierr); 7415582bec1SHong Zhang ierr = MatSetType(*newmat,MATSEQAIJ);CHKERRQ(ierr); 742e14861a4SHong Zhang 743e14861a4SHong Zhang nz_max = 0; 7445582bec1SHong Zhang for (i=0; i<m; i++) { 745e14861a4SHong Zhang nnz[i] = ml_cols[i+1] - ml_cols[i] + 1; 746e14861a4SHong Zhang if (nnz[i] > nz_max) nz_max = nnz[i]; 7475582bec1SHong Zhang } 7485582bec1SHong Zhang ierr = MatSeqAIJSetPreallocation(*newmat,0,nnz);CHKERRQ(ierr); 749*ab718edeSHong Zhang ierr = MatSetOption(*newmat,MAT_COLUMNS_SORTED);CHKERRQ(ierr); /* check! */ 7505582bec1SHong Zhang 751e14861a4SHong Zhang nz_max++; 752e14861a4SHong Zhang ierr = PetscMalloc(nz_max*(sizeof(int)+sizeof(PetscScalar)),&aj);CHKERRQ(ierr); 753e14861a4SHong Zhang aa = (PetscScalar*)(aj + nz_max); 75424a42b14SHong Zhang 7555582bec1SHong Zhang for (i=0; i<m; i++){ 756e14861a4SHong Zhang k = 0; 757e14861a4SHong Zhang /* diagonal entry */ 758e14861a4SHong Zhang aj[k] = i; aa[k++] = ml_vals[i]; 759e14861a4SHong Zhang /* off diagonal entries */ 760e14861a4SHong Zhang for (j=ml_cols[i]; j<ml_cols[i+1]; j++){ 761e14861a4SHong Zhang aj[k] = ml_cols[j]; aa[k++] = ml_vals[j]; 76224a42b14SHong Zhang } 763*ab718edeSHong Zhang /* sort aj and aa */ 764e14861a4SHong Zhang ierr = PetscSortIntWithScalarArray(nnz[i],aj,aa);CHKERRQ(ierr); 765e14861a4SHong Zhang ierr = MatSetValues(*newmat,1,&i,nnz[i],aj,aa,INSERT_VALUES);CHKERRQ(ierr); 7665582bec1SHong Zhang } 7675582bec1SHong Zhang ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 7685582bec1SHong Zhang ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 769e14861a4SHong Zhang ierr = PetscFree(aj);CHKERRQ(ierr); 7705582bec1SHong Zhang PetscFunctionReturn(0); 7715582bec1SHong Zhang } 7725582bec1SHong Zhang 7735582bec1SHong Zhang #undef __FUNCT__ 7745582bec1SHong Zhang #define __FUNCT__ "MatConvert_ML_SHELL" 7755582bec1SHong Zhang PetscErrorCode MatConvert_ML_SHELL(ML_Operator *mlmat,Mat *newmat) 7765582bec1SHong Zhang { 7775582bec1SHong Zhang PetscErrorCode ierr; 7785582bec1SHong Zhang PetscInt m,n; 7795582bec1SHong Zhang ML_Comm *MLcomm; 7805582bec1SHong Zhang Mat_MLShell *shellctx; 7815582bec1SHong Zhang 7825582bec1SHong Zhang PetscFunctionBegin; 7835582bec1SHong Zhang m = mlmat->outvec_leng; 7845582bec1SHong Zhang n = mlmat->invec_leng; 7855582bec1SHong Zhang if (!m || !n){ 7865582bec1SHong Zhang newmat = PETSC_NULL; 7875582bec1SHong Zhang } else { 7885582bec1SHong Zhang MLcomm = mlmat->comm; 7895582bec1SHong Zhang ierr = PetscNew(Mat_MLShell,&shellctx);CHKERRQ(ierr); 7905582bec1SHong Zhang ierr = PetscMemzero(shellctx,sizeof(Mat_MLShell));CHKERRQ(ierr); 7915582bec1SHong Zhang ierr = MatCreateShell(MLcomm->USR_comm,m,n,PETSC_DETERMINE,PETSC_DETERMINE,shellctx,newmat);CHKERRQ(ierr); 7925582bec1SHong Zhang ierr = MatShellSetOperation(*newmat,MATOP_MULT,(void *)MatMult_ML);CHKERRQ(ierr); 7935582bec1SHong Zhang ierr = MatShellSetOperation(*newmat,MATOP_MULT_ADD,(void *)MatMultAdd_ML);CHKERRQ(ierr); 7945582bec1SHong Zhang shellctx->A = *newmat; 7955582bec1SHong Zhang shellctx->mlmat = mlmat; 7965582bec1SHong Zhang ierr = VecCreate(PETSC_COMM_WORLD,&shellctx->y);CHKERRQ(ierr); 7975582bec1SHong Zhang ierr = VecSetSizes(shellctx->y,m,PETSC_DECIDE);CHKERRQ(ierr); 7985582bec1SHong Zhang ierr = VecSetFromOptions(shellctx->y);CHKERRQ(ierr); 7995582bec1SHong Zhang (*newmat)->ops->destroy = MatDestroy_ML; 8005582bec1SHong Zhang } 8015582bec1SHong Zhang PetscFunctionReturn(0); 8025582bec1SHong Zhang } 803e14861a4SHong Zhang 804e14861a4SHong Zhang #undef __FUNCT__ 805e14861a4SHong Zhang #define __FUNCT__ "MatConvert_ML_MPIAIJ" 806*ab718edeSHong Zhang PetscErrorCode MatConvert_ML_MPIAIJ(ML_Operator *mlmat,Mat *newmat) 807e14861a4SHong Zhang { 808*ab718edeSHong Zhang struct ML_CSR_MSRdata *matdata = (struct ML_CSR_MSRdata *)mlmat->data; 809*ab718edeSHong Zhang PetscInt *ml_cols=matdata->columns,*aj; 810*ab718edeSHong Zhang PetscScalar *ml_vals=matdata->values,*aa; 811e14861a4SHong Zhang PetscErrorCode ierr; 812*ab718edeSHong Zhang PetscInt i,j,k,*gordering; 813*ab718edeSHong Zhang PetscInt m=mlmat->outvec_leng,n,nnzA[m],nnzB[m],nnz[m],nz_max,row; 814*ab718edeSHong Zhang Mat A; 815e14861a4SHong Zhang 816e14861a4SHong Zhang PetscFunctionBegin; 817*ab718edeSHong Zhang if (mlmat->getrow == NULL) SETERRQ(PETSC_ERR_ARG_NULL,"mlmat->getrow = NULL"); 818*ab718edeSHong Zhang n = mlmat->invec_leng; 819e14861a4SHong Zhang if (m != n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"m %d must equal to n %d",m,n); 820e14861a4SHong Zhang 821*ab718edeSHong Zhang ierr = MatCreate(mlmat->comm->USR_comm,m,n,PETSC_DECIDE,PETSC_DECIDE,&A);CHKERRQ(ierr); 822*ab718edeSHong Zhang ierr = MatSetType(A,MATMPIAIJ);CHKERRQ(ierr); 823e14861a4SHong Zhang nz_max = 0; 824e14861a4SHong Zhang for (i=0; i<m; i++){ 825*ab718edeSHong Zhang nnz[i] = ml_cols[i+1] - ml_cols[i] + 1; 826e14861a4SHong Zhang if (nz_max < nnz[i]) nz_max = nnz[i]; 827*ab718edeSHong Zhang nnzA[i] = 1; /* diag */ 828*ab718edeSHong Zhang for (j=ml_cols[i]; j<ml_cols[i+1]; j++){ 829*ab718edeSHong Zhang if (ml_cols[j] < m) nnzA[i]++; 830e14861a4SHong Zhang } 831e14861a4SHong Zhang nnzB[i] = nnz[i] - nnzA[i]; 832e14861a4SHong Zhang } 833*ab718edeSHong Zhang ierr = MatMPIAIJSetPreallocation(A,0,nnzA,0,nnzB);CHKERRQ(ierr); 834e14861a4SHong Zhang 835*ab718edeSHong Zhang /* insert mat values -- remap row and column indices */ 836*ab718edeSHong Zhang nz_max++; 837*ab718edeSHong Zhang ierr = PetscMalloc(nz_max*(sizeof(int)+sizeof(PetscScalar)),&aj);CHKERRQ(ierr); 838*ab718edeSHong Zhang aa = (PetscScalar*)(aj + nz_max); 839*ab718edeSHong Zhang ML_build_global_numbering(mlmat,mlmat->comm,&gordering); 840e14861a4SHong Zhang for (i=0; i<m; i++){ 841e14861a4SHong Zhang row = gordering[i]; 842*ab718edeSHong Zhang k = 0; 843*ab718edeSHong Zhang /* diagonal entry */ 844*ab718edeSHong Zhang aj[k] = row; aa[k++] = ml_vals[i]; 845*ab718edeSHong Zhang /* off diagonal entries */ 846*ab718edeSHong Zhang for (j=ml_cols[i]; j<ml_cols[i+1]; j++){ 847*ab718edeSHong Zhang aj[k] = gordering[ml_cols[j]]; aa[k++] = ml_vals[j]; 848e14861a4SHong Zhang } 849*ab718edeSHong Zhang ierr = MatSetValues(A,1,&row,nnz[i],aj,aa,INSERT_VALUES);CHKERRQ(ierr); 850e14861a4SHong Zhang } 851*ab718edeSHong Zhang ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 852*ab718edeSHong Zhang ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 853*ab718edeSHong Zhang *newmat = A; 854e14861a4SHong Zhang 855*ab718edeSHong Zhang ierr = PetscFree(aj);CHKERRQ(ierr); 856e14861a4SHong Zhang PetscFunctionReturn(0); 857e14861a4SHong Zhang } 858