1dba47a55SKris Buschelman #define PETSCKSP_DLL 2ab718edeSHong Zhang 35582bec1SHong Zhang /* 41e5ab15bSHong Zhang Provides an interface to the ML 4.0 smoothed Aggregation 55582bec1SHong Zhang */ 66356e834SBarry Smith #include "private/pcimpl.h" /*I "petscpc.h" I*/ 75582bec1SHong Zhang #include "src/ksp/pc/impls/mg/mgimpl.h" /*I "petscmg.h" I*/ 85582bec1SHong Zhang #include "src/mat/impls/aij/seq/aij.h" 95582bec1SHong Zhang #include "src/mat/impls/aij/mpi/mpiaij.h" 10cb5d8e9eSHong Zhang 115582bec1SHong Zhang EXTERN_C_BEGIN 125582bec1SHong Zhang #include <math.h> 135582bec1SHong Zhang #include "ml_include.h" 145582bec1SHong Zhang EXTERN_C_END 155582bec1SHong Zhang 165582bec1SHong Zhang /* The context (data structure) at each grid level */ 175582bec1SHong Zhang typedef struct { 185582bec1SHong Zhang Vec x,b,r; /* global vectors */ 195582bec1SHong Zhang Mat A,P,R; 205582bec1SHong Zhang KSP ksp; 215582bec1SHong Zhang } GridCtx; 225582bec1SHong Zhang 235582bec1SHong Zhang /* The context used to input PETSc matrix into ML at fine grid */ 245582bec1SHong Zhang typedef struct { 255582bec1SHong Zhang Mat A,Aloc; 2624a42b14SHong Zhang Vec x,y; 275582bec1SHong Zhang ML_Operator *mlmat; 285582bec1SHong Zhang PetscScalar *pwork; /* tmp array used by PetscML_comm() */ 295582bec1SHong Zhang } FineGridCtx; 305582bec1SHong Zhang 315582bec1SHong Zhang /* The context associates a ML matrix with a PETSc shell matrix */ 325582bec1SHong Zhang typedef struct { 335582bec1SHong Zhang Mat A; /* PETSc shell matrix associated with mlmat */ 345582bec1SHong Zhang ML_Operator *mlmat; /* ML matrix assorciated with A */ 355582bec1SHong Zhang Vec y; 365582bec1SHong Zhang } Mat_MLShell; 375582bec1SHong Zhang 385582bec1SHong Zhang /* Private context for the ML preconditioner */ 395582bec1SHong Zhang typedef struct { 405582bec1SHong Zhang ML *ml_object; 415582bec1SHong Zhang ML_Aggregate *agg_object; 425582bec1SHong Zhang GridCtx *gridctx; 435582bec1SHong Zhang FineGridCtx *PetscMLdata; 445582bec1SHong Zhang PetscInt fine_level,MaxNlevels,MaxCoarseSize,CoarsenScheme; 455582bec1SHong Zhang PetscReal Threshold,DampingFactor; 465582bec1SHong Zhang PetscTruth SpectralNormScheme_Anorm; 475582bec1SHong Zhang PetscMPIInt size; 485582bec1SHong Zhang PetscErrorCode (*PCSetUp)(PC); 495582bec1SHong Zhang PetscErrorCode (*PCDestroy)(PC); 505582bec1SHong Zhang } PC_ML; 5141ca0015SHong Zhang 5241ca0015SHong Zhang extern int PetscML_getrow(ML_Operator *ML_data,int N_requested_rows,int requested_rows[], 535582bec1SHong Zhang int allocated_space,int columns[],double values[],int row_lengths[]); 5441ca0015SHong Zhang extern int PetscML_matvec(ML_Operator *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); 5875179d2cSHong Zhang extern PetscErrorCode MatConvert_MPIAIJ_ML(Mat,MatType,MatReuse,Mat*); 595582bec1SHong Zhang extern PetscErrorCode MatDestroy_ML(Mat); 60eef31507SHong Zhang extern PetscErrorCode MatWrapML_SeqAIJ(ML_Operator*,Mat*); 61eef31507SHong Zhang extern PetscErrorCode MatWrapML_MPIAIJ(ML_Operator*,Mat*); 62eef31507SHong Zhang extern PetscErrorCode MatWrapML_SHELL(ML_Operator*,Mat*); 635582bec1SHong Zhang 645582bec1SHong Zhang /* -------------------------------------------------------------------------- */ 655582bec1SHong Zhang /* 665582bec1SHong Zhang PCSetUp_ML - Prepares for the use of the ML preconditioner 675582bec1SHong Zhang by setting data structures and options. 685582bec1SHong Zhang 695582bec1SHong Zhang Input Parameter: 705582bec1SHong Zhang . pc - the preconditioner context 715582bec1SHong Zhang 725582bec1SHong Zhang Application Interface Routine: PCSetUp() 735582bec1SHong Zhang 745582bec1SHong Zhang Notes: 755582bec1SHong Zhang The interface routine PCSetUp() is not usually called directly by 765582bec1SHong Zhang the user, but instead is called by PCApply() if necessary. 775582bec1SHong Zhang */ 786ca4d86aSHong Zhang extern PetscErrorCode PCSetFromOptions_MG(PC); 795582bec1SHong Zhang #undef __FUNCT__ 805582bec1SHong Zhang #define __FUNCT__ "PCSetUp_ML" 816ca4d86aSHong Zhang PetscErrorCode PCSetUp_ML(PC pc) 825582bec1SHong Zhang { 835582bec1SHong Zhang PetscErrorCode ierr; 84eef31507SHong Zhang PetscMPIInt size; 855582bec1SHong Zhang FineGridCtx *PetscMLdata; 865582bec1SHong Zhang ML *ml_object; 875582bec1SHong Zhang ML_Aggregate *agg_object; 885582bec1SHong Zhang ML_Operator *mlmat; 89ac346b81SHong Zhang PetscInt nlocal_allcols,Nlevels,mllevel,level,level1,m,fine_level; 905582bec1SHong Zhang Mat A,Aloc; 915582bec1SHong Zhang GridCtx *gridctx; 925582bec1SHong Zhang PC_ML *pc_ml=PETSC_NULL; 93776b82aeSLisandro Dalcin PetscContainer container; 945582bec1SHong Zhang 955582bec1SHong Zhang PetscFunctionBegin; 965582bec1SHong Zhang ierr = PetscObjectQuery((PetscObject)pc,"PC_ML",(PetscObject *)&container);CHKERRQ(ierr); 975582bec1SHong Zhang if (container) { 98776b82aeSLisandro Dalcin ierr = PetscContainerGetPointer(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 pc_ml->size = size; 1095582bec1SHong Zhang if (size > 1){ 110eef31507SHong Zhang ierr = MatConvert_MPIAIJ_ML(A,PETSC_NULL,MAT_INITIAL_MATRIX,&Aloc);CHKERRQ(ierr); 1115582bec1SHong Zhang } else { 1125582bec1SHong Zhang Aloc = A; 1135582bec1SHong Zhang } 1145582bec1SHong Zhang 1155582bec1SHong Zhang /* create and initialize struct 'PetscMLdata' */ 1165582bec1SHong Zhang ierr = PetscNew(FineGridCtx,&PetscMLdata);CHKERRQ(ierr); 1175582bec1SHong Zhang PetscMLdata->A = A; 1185582bec1SHong Zhang PetscMLdata->Aloc = Aloc; 119a803a431SHong Zhang ierr = PetscMalloc((Aloc->cmap.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){ 124a803a431SHong Zhang ierr = VecSetSizes(PetscMLdata->x,A->cmap.n,A->cmap.n);CHKERRQ(ierr); 12524a42b14SHong Zhang } else { 126a803a431SHong Zhang ierr = VecSetSizes(PetscMLdata->x,Aloc->cmap.n,Aloc->cmap.n);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); 131a803a431SHong Zhang ierr = VecSetSizes(PetscMLdata->y,A->rmap.n,PETSC_DECIDE);CHKERRQ(ierr); 13224a42b14SHong Zhang ierr = VecSetType(PetscMLdata->y,VECSEQ);CHKERRQ(ierr); 13324a42b14SHong Zhang 1345582bec1SHong Zhang /* create ML discretization matrix at fine grid */ 135*45cf47abSHong Zhang /* ML requires input of fine-grid matrix. It determines nlevels. */ 1365582bec1SHong Zhang ierr = MatGetSize(Aloc,&m,&nlocal_allcols);CHKERRQ(ierr); 1375582bec1SHong Zhang ML_Create(&ml_object,pc_ml->MaxNlevels); 1385582bec1SHong Zhang ML_Init_Amatrix(ml_object,0,m,m,PetscMLdata); 1395582bec1SHong Zhang ML_Set_Amatrix_Getrow(ml_object,0,PetscML_getrow,PetscML_comm,nlocal_allcols); 1405582bec1SHong Zhang ML_Set_Amatrix_Matvec(ml_object,0,PetscML_matvec); 1415582bec1SHong Zhang 1425582bec1SHong Zhang /* aggregation */ 1435582bec1SHong Zhang ML_Aggregate_Create(&agg_object); 1445582bec1SHong Zhang ML_Aggregate_Set_MaxCoarseSize(agg_object,pc_ml->MaxCoarseSize); 1455582bec1SHong Zhang /* set options */ 1465582bec1SHong Zhang switch (pc_ml->CoarsenScheme) { 1475582bec1SHong Zhang case 1: 1485582bec1SHong Zhang ML_Aggregate_Set_CoarsenScheme_Coupled(agg_object);break; 1495582bec1SHong Zhang case 2: 1505582bec1SHong Zhang ML_Aggregate_Set_CoarsenScheme_MIS(agg_object);break; 1515582bec1SHong Zhang case 3: 1525582bec1SHong Zhang ML_Aggregate_Set_CoarsenScheme_METIS(agg_object);break; 1535582bec1SHong Zhang } 1545582bec1SHong Zhang ML_Aggregate_Set_Threshold(agg_object,pc_ml->Threshold); 1555582bec1SHong Zhang ML_Aggregate_Set_DampingFactor(agg_object,pc_ml->DampingFactor); 1565582bec1SHong Zhang if (pc_ml->SpectralNormScheme_Anorm){ 1575582bec1SHong Zhang ML_Aggregate_Set_SpectralNormScheme_Anorm(agg_object); 1585582bec1SHong Zhang } 1595582bec1SHong Zhang 1605582bec1SHong Zhang Nlevels = ML_Gen_MGHierarchy_UsingAggregation(ml_object,0,ML_INCREASING,agg_object); 1615582bec1SHong Zhang if (Nlevels<=0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Nlevels %d must > 0",Nlevels); 16297177400SBarry Smith ierr = PCMGSetLevels(pc,Nlevels,PETSC_NULL);CHKERRQ(ierr); 16397177400SBarry Smith ierr = PCSetFromOptions_MG(pc);CHKERRQ(ierr); /* should be called in PCSetFromOptions_ML(), but cannot be called prior to PCMGSetLevels() */ 1645582bec1SHong Zhang pc_ml->ml_object = ml_object; 1655582bec1SHong Zhang pc_ml->agg_object = agg_object; 1665582bec1SHong Zhang 1675582bec1SHong Zhang ierr = PetscMalloc(Nlevels*sizeof(GridCtx),&gridctx);CHKERRQ(ierr); 1685582bec1SHong Zhang fine_level = Nlevels - 1; 1695582bec1SHong Zhang pc_ml->gridctx = gridctx; 1705582bec1SHong Zhang pc_ml->fine_level = fine_level; 1715582bec1SHong Zhang 1725582bec1SHong Zhang /* wrap ML matrices by PETSc shell matrices at coarsened grids. 1735582bec1SHong Zhang Level 0 is the finest grid for ML, but coarsest for PETSc! */ 174e14861a4SHong Zhang gridctx[fine_level].A = A; 175e14861a4SHong Zhang level = fine_level - 1; 176ab718edeSHong Zhang if (size == 1){ /* convert ML P, R and A into seqaij format */ 1775582bec1SHong Zhang for (mllevel=1; mllevel<Nlevels; mllevel++){ 178e14861a4SHong Zhang mlmat = &(ml_object->Pmat[mllevel]); 179eef31507SHong Zhang ierr = MatWrapML_SeqAIJ(mlmat,&gridctx[level].P);CHKERRQ(ierr); 180e14861a4SHong Zhang mlmat = &(ml_object->Amat[mllevel]); 181eef31507SHong Zhang ierr = MatWrapML_SeqAIJ(mlmat,&gridctx[level].A);CHKERRQ(ierr); 182e14861a4SHong Zhang mlmat = &(ml_object->Rmat[mllevel-1]); 183eef31507SHong Zhang ierr = MatWrapML_SeqAIJ(mlmat,&gridctx[level].R);CHKERRQ(ierr); 1845582bec1SHong Zhang level--; 1855582bec1SHong Zhang } 186ab718edeSHong Zhang } else { /* convert ML P and R into shell format, ML A into mpiaij format */ 1875582bec1SHong Zhang for (mllevel=1; mllevel<Nlevels; mllevel++){ 1885582bec1SHong Zhang mlmat = &(ml_object->Pmat[mllevel]); 189eef31507SHong Zhang ierr = MatWrapML_SHELL(mlmat,&gridctx[level].P);CHKERRQ(ierr); 190ab718edeSHong Zhang mlmat = &(ml_object->Rmat[mllevel-1]); 191eef31507SHong Zhang ierr = MatWrapML_SHELL(mlmat,&gridctx[level].R);CHKERRQ(ierr); 1925582bec1SHong Zhang mlmat = &(ml_object->Amat[mllevel]); 193eef31507SHong Zhang ierr = MatWrapML_MPIAIJ(mlmat,&gridctx[level].A);CHKERRQ(ierr); 1945582bec1SHong Zhang level--; 1955582bec1SHong Zhang } 1965582bec1SHong Zhang } 1975582bec1SHong Zhang 1985582bec1SHong Zhang /* create coarse level and the interpolation between the levels */ 199ac346b81SHong Zhang for (level=0; level<fine_level; level++){ 2005582bec1SHong Zhang ierr = VecCreate(gridctx[level].A->comm,&gridctx[level].x);CHKERRQ(ierr); 201a803a431SHong Zhang ierr = VecSetSizes(gridctx[level].x,gridctx[level].A->cmap.n,PETSC_DECIDE);CHKERRQ(ierr); 2025582bec1SHong Zhang ierr = VecSetType(gridctx[level].x,VECMPI);CHKERRQ(ierr); 20397177400SBarry Smith ierr = PCMGSetX(pc,level,gridctx[level].x);CHKERRQ(ierr); 2045582bec1SHong Zhang 2055582bec1SHong Zhang ierr = VecCreate(gridctx[level].A->comm,&gridctx[level].b);CHKERRQ(ierr); 206a803a431SHong Zhang ierr = VecSetSizes(gridctx[level].b,gridctx[level].A->rmap.n,PETSC_DECIDE);CHKERRQ(ierr); 2075582bec1SHong Zhang ierr = VecSetType(gridctx[level].b,VECMPI);CHKERRQ(ierr); 20897177400SBarry Smith ierr = PCMGSetRhs(pc,level,gridctx[level].b);CHKERRQ(ierr); 209ac346b81SHong Zhang 210ac346b81SHong Zhang level1 = level + 1; 211ac346b81SHong Zhang ierr = VecCreate(gridctx[level1].A->comm,&gridctx[level1].r);CHKERRQ(ierr); 212a803a431SHong Zhang ierr = VecSetSizes(gridctx[level1].r,gridctx[level1].A->rmap.n,PETSC_DECIDE);CHKERRQ(ierr); 213ac346b81SHong Zhang ierr = VecSetType(gridctx[level1].r,VECMPI);CHKERRQ(ierr); 21497177400SBarry Smith ierr = PCMGSetR(pc,level1,gridctx[level1].r);CHKERRQ(ierr); 215ac346b81SHong Zhang 21697177400SBarry Smith ierr = PCMGSetInterpolate(pc,level1,gridctx[level].P);CHKERRQ(ierr); 21797177400SBarry Smith ierr = PCMGSetRestriction(pc,level1,gridctx[level].R);CHKERRQ(ierr); 2185582bec1SHong Zhang 2195582bec1SHong Zhang if (level == 0){ 22097177400SBarry Smith ierr = PCMGGetCoarseSolve(pc,&gridctx[level].ksp);CHKERRQ(ierr); 2215582bec1SHong Zhang } else { 22297177400SBarry Smith ierr = PCMGGetSmoother(pc,level,&gridctx[level].ksp);CHKERRQ(ierr); 22397177400SBarry Smith ierr = PCMGSetResidual(pc,level,PCMGDefaultResidual,gridctx[level].A);CHKERRQ(ierr); 2245582bec1SHong Zhang } 2255582bec1SHong Zhang ierr = KSPSetOperators(gridctx[level].ksp,gridctx[level].A,gridctx[level].A,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 2265582bec1SHong Zhang } 22797177400SBarry Smith ierr = PCMGGetSmoother(pc,fine_level,&gridctx[fine_level].ksp);CHKERRQ(ierr); 22897177400SBarry Smith ierr = PCMGSetResidual(pc,fine_level,PCMGDefaultResidual,gridctx[fine_level].A);CHKERRQ(ierr); 229ac346b81SHong Zhang ierr = KSPSetOperators(gridctx[fine_level].ksp,gridctx[level].A,gridctx[fine_level].A,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 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__ 238776b82aeSLisandro Dalcin #define __FUNCT__ "PetscContainerDestroy_PC_ML" 239776b82aeSLisandro Dalcin PetscErrorCode PetscContainerDestroy_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); 251ac346b81SHong Zhang if (pc_ml->PetscMLdata->x){ierr = VecDestroy(pc_ml->PetscMLdata->x);CHKERRQ(ierr);} 252ac346b81SHong Zhang if (pc_ml->PetscMLdata->y){ierr = VecDestroy(pc_ml->PetscMLdata->y);CHKERRQ(ierr);} 2535582bec1SHong Zhang ierr = PetscFree(pc_ml->PetscMLdata);CHKERRQ(ierr); 2545582bec1SHong Zhang 255ac346b81SHong Zhang for (level=0; level<pc_ml->fine_level; level++){ 256ac346b81SHong Zhang if (pc_ml->gridctx[level].A){ierr = MatDestroy(pc_ml->gridctx[level].A);CHKERRQ(ierr);} 257ac346b81SHong Zhang if (pc_ml->gridctx[level].P){ierr = MatDestroy(pc_ml->gridctx[level].P);CHKERRQ(ierr);} 258ac346b81SHong Zhang if (pc_ml->gridctx[level].R){ierr = MatDestroy(pc_ml->gridctx[level].R);CHKERRQ(ierr);} 259ac346b81SHong Zhang if (pc_ml->gridctx[level].x){ierr = VecDestroy(pc_ml->gridctx[level].x);CHKERRQ(ierr);} 260ac346b81SHong Zhang if (pc_ml->gridctx[level].b){ierr = VecDestroy(pc_ml->gridctx[level].b);CHKERRQ(ierr);} 261ac346b81SHong Zhang if (pc_ml->gridctx[level+1].r){ierr = VecDestroy(pc_ml->gridctx[level+1].r);CHKERRQ(ierr);} 2625582bec1SHong Zhang } 2635582bec1SHong Zhang ierr = PetscFree(pc_ml->gridctx);CHKERRQ(ierr); 2645582bec1SHong Zhang ierr = PetscFree(pc_ml);CHKERRQ(ierr); 2655582bec1SHong Zhang PetscFunctionReturn(0); 2665582bec1SHong Zhang } 2675582bec1SHong Zhang /* -------------------------------------------------------------------------- */ 2685582bec1SHong Zhang /* 2695582bec1SHong Zhang PCDestroy_ML - Destroys the private context for the ML preconditioner 2705582bec1SHong Zhang that was created with PCCreate_ML(). 2715582bec1SHong Zhang 2725582bec1SHong Zhang Input Parameter: 2735582bec1SHong Zhang . pc - the preconditioner context 2745582bec1SHong Zhang 2755582bec1SHong Zhang Application Interface Routine: PCDestroy() 2765582bec1SHong Zhang */ 2775582bec1SHong Zhang #undef __FUNCT__ 2785582bec1SHong Zhang #define __FUNCT__ "PCDestroy_ML" 2796ca4d86aSHong Zhang PetscErrorCode PCDestroy_ML(PC pc) 2805582bec1SHong Zhang { 2815582bec1SHong Zhang PetscErrorCode ierr; 2825582bec1SHong Zhang PC_ML *pc_ml=PETSC_NULL; 283776b82aeSLisandro Dalcin PetscContainer container; 2845582bec1SHong Zhang 2855582bec1SHong Zhang PetscFunctionBegin; 2865582bec1SHong Zhang ierr = PetscObjectQuery((PetscObject)pc,"PC_ML",(PetscObject *)&container);CHKERRQ(ierr); 2875582bec1SHong Zhang if (container) { 288776b82aeSLisandro Dalcin ierr = PetscContainerGetPointer(container,(void **)&pc_ml);CHKERRQ(ierr); 2895582bec1SHong Zhang pc->ops->destroy = pc_ml->PCDestroy; 2905582bec1SHong Zhang } else { 2915582bec1SHong Zhang SETERRQ(PETSC_ERR_ARG_NULL,"Container does not exit"); 2925582bec1SHong Zhang } 2939cb0ec93SHong Zhang /* detach pc and PC_ML and dereference container */ 2945582bec1SHong Zhang ierr = PetscObjectCompose((PetscObject)pc,"PC_ML",0);CHKERRQ(ierr); 2955582bec1SHong Zhang ierr = (*pc->ops->destroy)(pc);CHKERRQ(ierr); 2965582bec1SHong Zhang 297776b82aeSLisandro Dalcin ierr = PetscContainerDestroy(container);CHKERRQ(ierr); 2985582bec1SHong Zhang PetscFunctionReturn(0); 2995582bec1SHong Zhang } 3005582bec1SHong Zhang 3015582bec1SHong Zhang #undef __FUNCT__ 3025582bec1SHong Zhang #define __FUNCT__ "PCSetFromOptions_ML" 3036ca4d86aSHong Zhang PetscErrorCode PCSetFromOptions_ML(PC pc) 3045582bec1SHong Zhang { 3055582bec1SHong Zhang PetscErrorCode ierr; 3065582bec1SHong Zhang PetscInt indx,m,PrintLevel,MaxNlevels,MaxCoarseSize; 3075582bec1SHong Zhang PetscReal Threshold,DampingFactor; 3085582bec1SHong Zhang PetscTruth flg; 3095582bec1SHong Zhang const char *scheme[] = {"Uncoupled","Coupled","MIS","METIS"}; 3105582bec1SHong Zhang PC_ML *pc_ml=PETSC_NULL; 311776b82aeSLisandro Dalcin PetscContainer container; 3129dcbbd2bSBarry Smith PCMGType mgtype; 3135582bec1SHong Zhang 3145582bec1SHong Zhang PetscFunctionBegin; 3155582bec1SHong Zhang ierr = PetscObjectQuery((PetscObject)pc,"PC_ML",(PetscObject *)&container);CHKERRQ(ierr); 3165582bec1SHong Zhang if (container) { 317776b82aeSLisandro Dalcin ierr = PetscContainerGetPointer(container,(void **)&pc_ml);CHKERRQ(ierr); 3185582bec1SHong Zhang } else { 3195582bec1SHong Zhang SETERRQ(PETSC_ERR_ARG_NULL,"Container does not exit"); 3205582bec1SHong Zhang } 3216ca4d86aSHong Zhang 3225582bec1SHong Zhang /* inherited MG options */ 3236ca4d86aSHong Zhang ierr = PetscOptionsHead("Multigrid options(inherited)");CHKERRQ(ierr); 3245582bec1SHong Zhang ierr = PetscOptionsInt("-pc_mg_cycles","1 for V cycle, 2 for W-cycle","MGSetCycles",1,&m,&flg);CHKERRQ(ierr); 3255582bec1SHong Zhang ierr = PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","MGSetNumberSmoothUp",1,&m,&flg);CHKERRQ(ierr); 3265582bec1SHong Zhang ierr = PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","MGSetNumberSmoothDown",1,&m,&flg);CHKERRQ(ierr); 3279dcbbd2bSBarry Smith ierr = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)PC_MG_MULTIPLICATIVE,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr); 3285582bec1SHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 3295582bec1SHong Zhang 3305582bec1SHong Zhang /* ML options */ 3315582bec1SHong Zhang ierr = PetscOptionsHead("ML options");CHKERRQ(ierr); 3325582bec1SHong Zhang /* set defaults */ 3335582bec1SHong Zhang PrintLevel = 0; 3345582bec1SHong Zhang MaxNlevels = 10; 3355582bec1SHong Zhang MaxCoarseSize = 1; 3365582bec1SHong Zhang indx = 0; 3375582bec1SHong Zhang Threshold = 0.0; 3385582bec1SHong Zhang DampingFactor = 4.0/3.0; 3395582bec1SHong Zhang 3405582bec1SHong Zhang ierr = PetscOptionsInt("-pc_ml_PrintLevel","Print level","ML_Set_PrintLevel",PrintLevel,&PrintLevel,PETSC_NULL);CHKERRQ(ierr); 3415582bec1SHong Zhang ML_Set_PrintLevel(PrintLevel); 3425582bec1SHong Zhang 3435582bec1SHong Zhang ierr = PetscOptionsInt("-pc_ml_maxNlevels","Maximum number of levels","None",MaxNlevels,&MaxNlevels,PETSC_NULL);CHKERRQ(ierr); 3445582bec1SHong Zhang pc_ml->MaxNlevels = MaxNlevels; 3455582bec1SHong Zhang 3465582bec1SHong Zhang ierr = PetscOptionsInt("-pc_ml_maxCoarseSize","Maximum coarsest mesh size","ML_Aggregate_Set_MaxCoarseSize",MaxCoarseSize,&MaxCoarseSize,PETSC_NULL);CHKERRQ(ierr); 3475582bec1SHong Zhang pc_ml->MaxCoarseSize = MaxCoarseSize; 3485582bec1SHong Zhang 3495582bec1SHong Zhang ierr = PetscOptionsEList("-pc_ml_CoarsenScheme","Aggregate Coarsen Scheme","ML_Aggregate_Set_CoarsenScheme_*",scheme,4,scheme[0],&indx,PETSC_NULL);CHKERRQ(ierr); 3505582bec1SHong Zhang pc_ml->CoarsenScheme = indx; 3515582bec1SHong Zhang 3525582bec1SHong Zhang ierr = PetscOptionsReal("-pc_ml_DampingFactor","P damping factor","ML_Aggregate_Set_DampingFactor",DampingFactor,&DampingFactor,PETSC_NULL);CHKERRQ(ierr); 3535582bec1SHong Zhang pc_ml->DampingFactor = DampingFactor; 3545582bec1SHong Zhang 3555582bec1SHong Zhang ierr = PetscOptionsReal("-pc_ml_Threshold","Smoother drop tol","ML_Aggregate_Set_Threshold",Threshold,&Threshold,PETSC_NULL);CHKERRQ(ierr); 3565582bec1SHong Zhang pc_ml->Threshold = Threshold; 3575582bec1SHong Zhang 358a146b4dcSHong Zhang ierr = PetscOptionsTruth("-pc_ml_SpectralNormScheme_Anorm","Method used for estimating spectral radius","ML_Aggregate_Set_SpectralNormScheme_Anorm",PETSC_FALSE,&pc_ml->SpectralNormScheme_Anorm,PETSC_NULL); 3595582bec1SHong Zhang 3605582bec1SHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 3615582bec1SHong Zhang PetscFunctionReturn(0); 3625582bec1SHong Zhang } 3635582bec1SHong Zhang 3645582bec1SHong Zhang /* -------------------------------------------------------------------------- */ 3655582bec1SHong Zhang /* 3665582bec1SHong Zhang PCCreate_ML - Creates a ML preconditioner context, PC_ML, 3675582bec1SHong Zhang and sets this as the private data within the generic preconditioning 3685582bec1SHong Zhang context, PC, that was created within PCCreate(). 3695582bec1SHong Zhang 3705582bec1SHong Zhang Input Parameter: 3715582bec1SHong Zhang . pc - the preconditioner context 3725582bec1SHong Zhang 3735582bec1SHong Zhang Application Interface Routine: PCCreate() 3745582bec1SHong Zhang */ 3755582bec1SHong Zhang 3765582bec1SHong Zhang /*MC 3771e5ab15bSHong Zhang PCML - Use algebraic multigrid preconditioning. This preconditioner requires you provide 3785582bec1SHong Zhang fine grid discretization matrix. The coarser grid matrices and restriction/interpolation 3796ca4d86aSHong Zhang operators are computed by ML, with the matrices coverted to PETSc matrices in aij format 3806ca4d86aSHong Zhang and the restriction/interpolation operators wrapped as PETSc shell matrices. 3815582bec1SHong Zhang 3826ca4d86aSHong Zhang Options Database Key: 3836ca4d86aSHong Zhang Multigrid options(inherited) 3846ca4d86aSHong Zhang + -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (MGSetCycles) 3856ca4d86aSHong Zhang . -pc_mg_smoothup <1>: Number of post-smoothing steps (MGSetNumberSmoothUp) 3866ca4d86aSHong Zhang . -pc_mg_smoothdown <1>: Number of pre-smoothing steps (MGSetNumberSmoothDown) 3876ca4d86aSHong Zhang - -pc_mg_type <multiplicative> (one of) additive multiplicative full cascade kascade 3886ca4d86aSHong Zhang 38951f519a2SBarry Smith ML options: 3906ca4d86aSHong Zhang + -pc_ml_PrintLevel <0>: Print level (ML_Set_PrintLevel) 3916ca4d86aSHong Zhang . -pc_ml_maxNlevels <10>: Maximum number of levels (None) 3926ca4d86aSHong Zhang . -pc_ml_maxCoarseSize <1>: Maximum coarsest mesh size (ML_Aggregate_Set_MaxCoarseSize) 3936ca4d86aSHong Zhang . -pc_ml_CoarsenScheme <Uncoupled> (one of) Uncoupled Coupled MIS METIS 3946ca4d86aSHong Zhang . -pc_ml_DampingFactor <1.33333>: P damping factor (ML_Aggregate_Set_DampingFactor) 3956ca4d86aSHong Zhang . -pc_ml_Threshold <0>: Smoother drop tol (ML_Aggregate_Set_Threshold) 3966ca4d86aSHong Zhang - -pc_ml_SpectralNormScheme_Anorm: <false> Method used for estimating spectral radius (ML_Aggregate_Set_SpectralNormScheme_Anorm) 3975582bec1SHong Zhang 3985582bec1SHong Zhang Level: intermediate 3995582bec1SHong Zhang 4005582bec1SHong Zhang Concepts: multigrid 4015582bec1SHong Zhang 4025582bec1SHong Zhang .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, 40397177400SBarry Smith PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), MPSetCycles(), PCMGSetNumberSmoothDown(), 40497177400SBarry Smith PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(), 40597177400SBarry Smith PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(), 40697177400SBarry Smith PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR() 4075582bec1SHong Zhang M*/ 4085582bec1SHong Zhang 4095582bec1SHong Zhang EXTERN_C_BEGIN 4105582bec1SHong Zhang #undef __FUNCT__ 4115582bec1SHong Zhang #define __FUNCT__ "PCCreate_ML" 412dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_ML(PC pc) 4135582bec1SHong Zhang { 4145582bec1SHong Zhang PetscErrorCode ierr; 4155582bec1SHong Zhang PC_ML *pc_ml; 416776b82aeSLisandro Dalcin PetscContainer container; 4175582bec1SHong Zhang 4185582bec1SHong Zhang PetscFunctionBegin; 4195582bec1SHong Zhang /* initialize pc as PCMG */ 4205582bec1SHong Zhang ierr = PCSetType(pc,PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */ 4215582bec1SHong Zhang 4225582bec1SHong Zhang /* create a supporting struct and attach it to pc */ 4235582bec1SHong Zhang ierr = PetscNew(PC_ML,&pc_ml);CHKERRQ(ierr); 424776b82aeSLisandro Dalcin ierr = PetscContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr); 425776b82aeSLisandro Dalcin ierr = PetscContainerSetPointer(container,pc_ml);CHKERRQ(ierr); 426776b82aeSLisandro Dalcin ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_PC_ML);CHKERRQ(ierr); 4275582bec1SHong Zhang ierr = PetscObjectCompose((PetscObject)pc,"PC_ML",(PetscObject)container);CHKERRQ(ierr); 4285582bec1SHong Zhang 4295582bec1SHong Zhang pc_ml->PCSetUp = pc->ops->setup; 4305582bec1SHong Zhang pc_ml->PCDestroy = pc->ops->destroy; 4315582bec1SHong Zhang 4325582bec1SHong Zhang /* overwrite the pointers of PCMG by the functions of PCML */ 4335582bec1SHong Zhang pc->ops->setfromoptions = PCSetFromOptions_ML; 4345582bec1SHong Zhang pc->ops->setup = PCSetUp_ML; 4355582bec1SHong Zhang pc->ops->destroy = PCDestroy_ML; 4365582bec1SHong Zhang PetscFunctionReturn(0); 4375582bec1SHong Zhang } 4385582bec1SHong Zhang EXTERN_C_END 4395582bec1SHong Zhang 44041ca0015SHong Zhang int PetscML_getrow(ML_Operator *ML_data, int N_requested_rows, int requested_rows[], 4415582bec1SHong Zhang int allocated_space, int columns[], double values[], int row_lengths[]) 4425582bec1SHong Zhang { 4435582bec1SHong Zhang PetscErrorCode ierr; 4445582bec1SHong Zhang Mat Aloc; 4455582bec1SHong Zhang Mat_SeqAIJ *a; 4465582bec1SHong Zhang PetscInt m,i,j,k=0,row,*aj; 4475582bec1SHong Zhang PetscScalar *aa; 44841ca0015SHong Zhang FineGridCtx *ml=(FineGridCtx*)ML_Get_MyGetrowData(ML_data); 4495582bec1SHong Zhang 4505582bec1SHong Zhang Aloc = ml->Aloc; 4515582bec1SHong Zhang a = (Mat_SeqAIJ*)Aloc->data; 4525582bec1SHong Zhang ierr = MatGetSize(Aloc,&m,PETSC_NULL);CHKERRQ(ierr); 4535582bec1SHong Zhang 4545582bec1SHong Zhang for (i = 0; i<N_requested_rows; i++) { 4555582bec1SHong Zhang row = requested_rows[i]; 4565582bec1SHong Zhang row_lengths[i] = a->ilen[row]; 4575582bec1SHong Zhang if (allocated_space < k+row_lengths[i]) return(0); 4585582bec1SHong Zhang if ( (row >= 0) || (row <= (m-1)) ) { 4595582bec1SHong Zhang aj = a->j + a->i[row]; 4605582bec1SHong Zhang aa = a->a + a->i[row]; 4615582bec1SHong Zhang for (j=0; j<row_lengths[i]; j++){ 4625582bec1SHong Zhang columns[k] = aj[j]; 4635582bec1SHong Zhang values[k++] = aa[j]; 4645582bec1SHong Zhang } 4655582bec1SHong Zhang } 4665582bec1SHong Zhang } 4675582bec1SHong Zhang return(1); 4685582bec1SHong Zhang } 4695582bec1SHong Zhang 47041ca0015SHong Zhang int PetscML_matvec(ML_Operator *ML_data,int in_length,double p[],int out_length,double ap[]) 4715582bec1SHong Zhang { 4725582bec1SHong Zhang PetscErrorCode ierr; 47341ca0015SHong Zhang FineGridCtx *ml=(FineGridCtx*)ML_Get_MyMatvecData(ML_data); 4745582bec1SHong Zhang Mat A=ml->A, Aloc=ml->Aloc; 4755582bec1SHong Zhang PetscMPIInt size; 4765582bec1SHong Zhang PetscScalar *pwork=ml->pwork; 4775582bec1SHong Zhang PetscInt i; 4785582bec1SHong Zhang 4795582bec1SHong Zhang ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr); 4805582bec1SHong Zhang if (size == 1){ 48124a42b14SHong Zhang ierr = VecPlaceArray(ml->x,p);CHKERRQ(ierr); 4825582bec1SHong Zhang } else { 4835582bec1SHong Zhang for (i=0; i<in_length; i++) pwork[i] = p[i]; 4845582bec1SHong Zhang PetscML_comm(pwork,ml); 48524a42b14SHong Zhang ierr = VecPlaceArray(ml->x,pwork);CHKERRQ(ierr); 4865582bec1SHong Zhang } 48724a42b14SHong Zhang ierr = VecPlaceArray(ml->y,ap);CHKERRQ(ierr); 48824a42b14SHong Zhang ierr = MatMult(Aloc,ml->x,ml->y);CHKERRQ(ierr); 489957c941cSHong Zhang ierr = VecResetArray(ml->x);CHKERRQ(ierr); 490957c941cSHong Zhang ierr = VecResetArray(ml->y);CHKERRQ(ierr); 4915582bec1SHong Zhang return 0; 4925582bec1SHong Zhang } 4935582bec1SHong Zhang 4945582bec1SHong Zhang int PetscML_comm(double p[],void *ML_data) 4955582bec1SHong Zhang { 4965582bec1SHong Zhang PetscErrorCode ierr; 4975582bec1SHong Zhang FineGridCtx *ml=(FineGridCtx*)ML_data; 4985582bec1SHong Zhang Mat A=ml->A; 4995582bec1SHong Zhang Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data; 5005582bec1SHong Zhang PetscMPIInt size; 501a803a431SHong Zhang PetscInt i,in_length=A->rmap.n,out_length=ml->Aloc->cmap.n; 5025582bec1SHong Zhang PetscScalar *array; 5035582bec1SHong Zhang 5045582bec1SHong Zhang ierr = MPI_Comm_size(A->comm,&size);CHKERRQ(ierr); 5055582bec1SHong Zhang if (size == 1) return 0; 50624a42b14SHong Zhang 50724a42b14SHong Zhang ierr = VecPlaceArray(ml->y,p);CHKERRQ(ierr); 50824a42b14SHong Zhang ierr = VecScatterBegin(ml->y,a->lvec,INSERT_VALUES,SCATTER_FORWARD,a->Mvctx);CHKERRQ(ierr); 50924a42b14SHong Zhang ierr = VecScatterEnd(ml->y,a->lvec,INSERT_VALUES,SCATTER_FORWARD,a->Mvctx);CHKERRQ(ierr); 5107d15518fSHong Zhang ierr = VecResetArray(ml->y);CHKERRQ(ierr); 5115582bec1SHong Zhang ierr = VecGetArray(a->lvec,&array);CHKERRQ(ierr); 5125582bec1SHong Zhang for (i=in_length; i<out_length; i++){ 5135582bec1SHong Zhang p[i] = array[i-in_length]; 5145582bec1SHong Zhang } 5157d15518fSHong Zhang ierr = VecRestoreArray(a->lvec,&array);CHKERRQ(ierr); 5165582bec1SHong Zhang return 0; 5175582bec1SHong Zhang } 5185582bec1SHong Zhang #undef __FUNCT__ 5195582bec1SHong Zhang #define __FUNCT__ "MatMult_ML" 5205582bec1SHong Zhang PetscErrorCode MatMult_ML(Mat A,Vec x,Vec y) 5215582bec1SHong Zhang { 5225582bec1SHong Zhang PetscErrorCode ierr; 5235582bec1SHong Zhang Mat_MLShell *shell; 5245582bec1SHong Zhang PetscScalar *xarray,*yarray; 5255582bec1SHong Zhang PetscInt x_length,y_length; 5265582bec1SHong Zhang 5275582bec1SHong Zhang PetscFunctionBegin; 528a146b4dcSHong Zhang ierr = MatShellGetContext(A,(void **)&shell);CHKERRQ(ierr); 5295582bec1SHong Zhang ierr = VecGetArray(x,&xarray);CHKERRQ(ierr); 5305582bec1SHong Zhang ierr = VecGetArray(y,&yarray);CHKERRQ(ierr); 5315582bec1SHong Zhang x_length = shell->mlmat->invec_leng; 5325582bec1SHong Zhang y_length = shell->mlmat->outvec_leng; 5335582bec1SHong Zhang 5345582bec1SHong Zhang ML_Operator_Apply(shell->mlmat,x_length,xarray,y_length,yarray); 5355582bec1SHong Zhang 5365582bec1SHong Zhang ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr); 5375582bec1SHong Zhang ierr = VecRestoreArray(y,&yarray);CHKERRQ(ierr); 5385582bec1SHong Zhang PetscFunctionReturn(0); 5395582bec1SHong Zhang } 5405582bec1SHong Zhang /* MatMultAdd_ML - Compute y = w + A*x */ 5415582bec1SHong Zhang #undef __FUNCT__ 5425582bec1SHong Zhang #define __FUNCT__ "MatMultAdd_ML" 5435582bec1SHong Zhang PetscErrorCode MatMultAdd_ML(Mat A,Vec x,Vec w,Vec y) 5445582bec1SHong Zhang { 5455582bec1SHong Zhang PetscErrorCode ierr; 5465582bec1SHong Zhang Mat_MLShell *shell; 5475582bec1SHong Zhang PetscScalar *xarray,*yarray; 5485582bec1SHong Zhang PetscInt x_length,y_length; 5495582bec1SHong Zhang 5505582bec1SHong Zhang PetscFunctionBegin; 551a146b4dcSHong Zhang ierr = MatShellGetContext(A,(void **)&shell);CHKERRQ(ierr); 5525582bec1SHong Zhang ierr = VecGetArray(x,&xarray);CHKERRQ(ierr); 5535582bec1SHong Zhang ierr = VecGetArray(y,&yarray);CHKERRQ(ierr); 5545582bec1SHong Zhang 5555582bec1SHong Zhang x_length = shell->mlmat->invec_leng; 5565582bec1SHong Zhang y_length = shell->mlmat->outvec_leng; 5575582bec1SHong Zhang 5585582bec1SHong Zhang ML_Operator_Apply(shell->mlmat,x_length,xarray,y_length,yarray); 5595582bec1SHong Zhang 5605582bec1SHong Zhang ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr); 5615582bec1SHong Zhang ierr = VecRestoreArray(y,&yarray);CHKERRQ(ierr); 562efb30889SBarry Smith ierr = VecAXPY(y,1.0,w);CHKERRQ(ierr); 5635582bec1SHong Zhang 5645582bec1SHong Zhang PetscFunctionReturn(0); 5655582bec1SHong Zhang } 5665582bec1SHong Zhang 5675582bec1SHong Zhang /* newtype is ignored because "ml" is not listed under Petsc MatType yet */ 5685582bec1SHong Zhang #undef __FUNCT__ 5695582bec1SHong Zhang #define __FUNCT__ "MatConvert_MPIAIJ_ML" 57075179d2cSHong Zhang PetscErrorCode MatConvert_MPIAIJ_ML(Mat A,MatType newtype,MatReuse scall,Mat *Aloc) 5715582bec1SHong Zhang { 5725582bec1SHong Zhang PetscErrorCode ierr; 5735582bec1SHong Zhang Mat_MPIAIJ *mpimat=(Mat_MPIAIJ*)A->data; 5745582bec1SHong Zhang Mat_SeqAIJ *mat,*a=(Mat_SeqAIJ*)(mpimat->A)->data,*b=(Mat_SeqAIJ*)(mpimat->B)->data; 5755582bec1SHong Zhang PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j; 5765582bec1SHong Zhang PetscScalar *aa=a->a,*ba=b->a,*ca; 577a803a431SHong Zhang PetscInt am=A->rmap.n,an=A->cmap.n,i,j,k; 5785582bec1SHong Zhang PetscInt *ci,*cj,ncols; 5795582bec1SHong Zhang 5805582bec1SHong Zhang PetscFunctionBegin; 5815582bec1SHong Zhang if (am != an) SETERRQ2(PETSC_ERR_ARG_WRONG,"A must have a square diagonal portion, am: %d != an: %d",am,an); 5825582bec1SHong Zhang 5835582bec1SHong Zhang if (scall == MAT_INITIAL_MATRIX){ 5845582bec1SHong Zhang ierr = PetscMalloc((1+am)*sizeof(PetscInt),&ci);CHKERRQ(ierr); 5855582bec1SHong Zhang ci[0] = 0; 5865582bec1SHong Zhang for (i=0; i<am; i++){ 5875582bec1SHong Zhang ci[i+1] = ci[i] + (ai[i+1] - ai[i]) + (bi[i+1] - bi[i]); 5885582bec1SHong Zhang } 5895582bec1SHong Zhang ierr = PetscMalloc((1+ci[am])*sizeof(PetscInt),&cj);CHKERRQ(ierr); 5905582bec1SHong Zhang ierr = PetscMalloc((1+ci[am])*sizeof(PetscScalar),&ca);CHKERRQ(ierr); 5915582bec1SHong Zhang 5925582bec1SHong Zhang k = 0; 5935582bec1SHong Zhang for (i=0; i<am; i++){ 5945582bec1SHong Zhang /* diagonal portion of A */ 5955582bec1SHong Zhang ncols = ai[i+1] - ai[i]; 5965582bec1SHong Zhang for (j=0; j<ncols; j++) { 5975582bec1SHong Zhang cj[k] = *aj++; 5985582bec1SHong Zhang ca[k++] = *aa++; 5995582bec1SHong Zhang } 6005582bec1SHong Zhang /* off-diagonal portion of A */ 6015582bec1SHong Zhang ncols = bi[i+1] - bi[i]; 6025582bec1SHong Zhang for (j=0; j<ncols; j++) { 6035582bec1SHong Zhang cj[k] = an + (*bj); bj++; 6045582bec1SHong Zhang ca[k++] = *ba++; 6055582bec1SHong Zhang } 6065582bec1SHong Zhang } 6075582bec1SHong Zhang if (k != ci[am]) SETERRQ2(PETSC_ERR_ARG_WRONG,"k: %d != ci[am]: %d",k,ci[am]); 6085582bec1SHong Zhang 6095582bec1SHong Zhang /* put together the new matrix */ 610a803a431SHong Zhang an = mpimat->A->cmap.n+mpimat->B->cmap.n; 6115582bec1SHong Zhang ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,an,ci,cj,ca,Aloc);CHKERRQ(ierr); 6125582bec1SHong Zhang 6135582bec1SHong Zhang /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 6145582bec1SHong Zhang /* Since these are PETSc arrays, change flags to free them as necessary. */ 6155582bec1SHong Zhang mat = (Mat_SeqAIJ*)(*Aloc)->data; 6163756052fSSatish Balay mat->free_a = PETSC_TRUE; 6173756052fSSatish Balay mat->free_ij = PETSC_TRUE; 6183756052fSSatish Balay 6195582bec1SHong Zhang mat->nonew = 0; 6205582bec1SHong Zhang } else if (scall == MAT_REUSE_MATRIX){ 6215582bec1SHong Zhang mat=(Mat_SeqAIJ*)(*Aloc)->data; 6225582bec1SHong Zhang ci = mat->i; cj = mat->j; ca = mat->a; 6235582bec1SHong Zhang for (i=0; i<am; i++) { 6245582bec1SHong Zhang /* diagonal portion of A */ 6255582bec1SHong Zhang ncols = ai[i+1] - ai[i]; 6265582bec1SHong Zhang for (j=0; j<ncols; j++) *ca++ = *aa++; 6275582bec1SHong Zhang /* off-diagonal portion of A */ 6285582bec1SHong Zhang ncols = bi[i+1] - bi[i]; 6295582bec1SHong Zhang for (j=0; j<ncols; j++) *ca++ = *ba++; 6305582bec1SHong Zhang } 6315582bec1SHong Zhang } else { 6325582bec1SHong Zhang SETERRQ1(PETSC_ERR_ARG_WRONG,"Invalid MatReuse %d",(int)scall); 6335582bec1SHong Zhang } 6345582bec1SHong Zhang PetscFunctionReturn(0); 6355582bec1SHong Zhang } 6365582bec1SHong Zhang extern PetscErrorCode MatDestroy_Shell(Mat); 6375582bec1SHong Zhang #undef __FUNCT__ 6385582bec1SHong Zhang #define __FUNCT__ "MatDestroy_ML" 6395582bec1SHong Zhang PetscErrorCode MatDestroy_ML(Mat A) 6405582bec1SHong Zhang { 6415582bec1SHong Zhang PetscErrorCode ierr; 6425582bec1SHong Zhang Mat_MLShell *shell; 6435582bec1SHong Zhang 6445582bec1SHong Zhang PetscFunctionBegin; 645a146b4dcSHong Zhang ierr = MatShellGetContext(A,(void **)&shell);CHKERRQ(ierr); 6465582bec1SHong Zhang ierr = VecDestroy(shell->y);CHKERRQ(ierr); 6475582bec1SHong Zhang ierr = PetscFree(shell);CHKERRQ(ierr); 6485582bec1SHong Zhang ierr = MatDestroy_Shell(A);CHKERRQ(ierr); 649dbd8c25aSHong Zhang ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr); 6505582bec1SHong Zhang PetscFunctionReturn(0); 6515582bec1SHong Zhang } 6525582bec1SHong Zhang 6535582bec1SHong Zhang #undef __FUNCT__ 654eef31507SHong Zhang #define __FUNCT__ "MatWrapML_SeqAIJ" 655eef31507SHong Zhang PetscErrorCode MatWrapML_SeqAIJ(ML_Operator *mlmat,Mat *newmat) 6565582bec1SHong Zhang { 657e14861a4SHong Zhang struct ML_CSR_MSRdata *matdata = (struct ML_CSR_MSRdata *)mlmat->data; 6585582bec1SHong Zhang PetscErrorCode ierr; 6593e826d7bSSatish Balay PetscInt m=mlmat->outvec_leng,n=mlmat->invec_leng,*nnz,nz_max; 660e14861a4SHong Zhang PetscInt *ml_cols=matdata->columns,*aj,i,j,k; 661e14861a4SHong Zhang PetscScalar *ml_vals=matdata->values,*aa; 6625582bec1SHong Zhang 6635582bec1SHong Zhang PetscFunctionBegin; 664e14861a4SHong Zhang if ( mlmat->getrow == NULL) SETERRQ(PETSC_ERR_ARG_NULL,"mlmat->getrow = NULL"); 665ab718edeSHong Zhang if (m != n){ /* ML Pmat and Rmat are in CSR format. Pass array pointers into SeqAIJ matrix */ 666e14861a4SHong Zhang ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,m,n,matdata->rowptr,ml_cols,ml_vals,newmat);CHKERRQ(ierr); 66724a42b14SHong Zhang PetscFunctionReturn(0); 66824a42b14SHong Zhang } 6695582bec1SHong Zhang 670e14861a4SHong Zhang /* ML Amat is in MSR format. Copy its data into SeqAIJ matrix */ 671f69a0ea3SMatthew Knepley ierr = MatCreate(PETSC_COMM_SELF,newmat);CHKERRQ(ierr); 672f69a0ea3SMatthew Knepley ierr = MatSetSizes(*newmat,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 6735582bec1SHong Zhang ierr = MatSetType(*newmat,MATSEQAIJ);CHKERRQ(ierr); 6743e826d7bSSatish Balay ierr = PetscMalloc((m+1)*sizeof(PetscInt),&nnz); 675e14861a4SHong Zhang 676e14861a4SHong Zhang nz_max = 0; 6775582bec1SHong Zhang for (i=0; i<m; i++) { 678e14861a4SHong Zhang nnz[i] = ml_cols[i+1] - ml_cols[i] + 1; 679e14861a4SHong Zhang if (nnz[i] > nz_max) nz_max = nnz[i]; 6805582bec1SHong Zhang } 6815582bec1SHong Zhang ierr = MatSeqAIJSetPreallocation(*newmat,0,nnz);CHKERRQ(ierr); 682ab718edeSHong Zhang ierr = MatSetOption(*newmat,MAT_COLUMNS_SORTED);CHKERRQ(ierr); /* check! */ 6835582bec1SHong Zhang 684e14861a4SHong Zhang nz_max++; 6857c307921SBarry Smith ierr = PetscMalloc(nz_max*(sizeof(PetscInt)+sizeof(PetscScalar)),&aj);CHKERRQ(ierr); 686e14861a4SHong Zhang aa = (PetscScalar*)(aj + nz_max); 68724a42b14SHong Zhang 6885582bec1SHong Zhang for (i=0; i<m; i++){ 689e14861a4SHong Zhang k = 0; 690e14861a4SHong Zhang /* diagonal entry */ 691e14861a4SHong Zhang aj[k] = i; aa[k++] = ml_vals[i]; 692e14861a4SHong Zhang /* off diagonal entries */ 693e14861a4SHong Zhang for (j=ml_cols[i]; j<ml_cols[i+1]; j++){ 694e14861a4SHong Zhang aj[k] = ml_cols[j]; aa[k++] = ml_vals[j]; 69524a42b14SHong Zhang } 696ab718edeSHong Zhang /* sort aj and aa */ 697e14861a4SHong Zhang ierr = PetscSortIntWithScalarArray(nnz[i],aj,aa);CHKERRQ(ierr); 698e14861a4SHong Zhang ierr = MatSetValues(*newmat,1,&i,nnz[i],aj,aa,INSERT_VALUES);CHKERRQ(ierr); 6995582bec1SHong Zhang } 7005582bec1SHong Zhang ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 7015582bec1SHong Zhang ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 702e14861a4SHong Zhang ierr = PetscFree(aj);CHKERRQ(ierr); 7033e826d7bSSatish Balay ierr = PetscFree(nnz);CHKERRQ(ierr); 7045582bec1SHong Zhang PetscFunctionReturn(0); 7055582bec1SHong Zhang } 7065582bec1SHong Zhang 7075582bec1SHong Zhang #undef __FUNCT__ 708eef31507SHong Zhang #define __FUNCT__ "MatWrapML_SHELL" 709eef31507SHong Zhang PetscErrorCode MatWrapML_SHELL(ML_Operator *mlmat,Mat *newmat) 7105582bec1SHong Zhang { 7115582bec1SHong Zhang PetscErrorCode ierr; 7125582bec1SHong Zhang PetscInt m,n; 7135582bec1SHong Zhang ML_Comm *MLcomm; 7145582bec1SHong Zhang Mat_MLShell *shellctx; 7155582bec1SHong Zhang 7165582bec1SHong Zhang PetscFunctionBegin; 7175582bec1SHong Zhang m = mlmat->outvec_leng; 7185582bec1SHong Zhang n = mlmat->invec_leng; 7195582bec1SHong Zhang if (!m || !n){ 7205582bec1SHong Zhang newmat = PETSC_NULL; 7215582bec1SHong Zhang } else { 7225582bec1SHong Zhang MLcomm = mlmat->comm; 7235582bec1SHong Zhang ierr = PetscNew(Mat_MLShell,&shellctx);CHKERRQ(ierr); 7245582bec1SHong Zhang ierr = MatCreateShell(MLcomm->USR_comm,m,n,PETSC_DETERMINE,PETSC_DETERMINE,shellctx,newmat);CHKERRQ(ierr); 7253e826d7bSSatish Balay ierr = MatShellSetOperation(*newmat,MATOP_MULT,(void(*)(void))MatMult_ML);CHKERRQ(ierr); 7263e826d7bSSatish Balay ierr = MatShellSetOperation(*newmat,MATOP_MULT_ADD,(void(*)(void))MatMultAdd_ML);CHKERRQ(ierr); 7275582bec1SHong Zhang shellctx->A = *newmat; 7285582bec1SHong Zhang shellctx->mlmat = mlmat; 7295582bec1SHong Zhang ierr = VecCreate(PETSC_COMM_WORLD,&shellctx->y);CHKERRQ(ierr); 7305582bec1SHong Zhang ierr = VecSetSizes(shellctx->y,m,PETSC_DECIDE);CHKERRQ(ierr); 7315582bec1SHong Zhang ierr = VecSetFromOptions(shellctx->y);CHKERRQ(ierr); 7325582bec1SHong Zhang (*newmat)->ops->destroy = MatDestroy_ML; 7335582bec1SHong Zhang } 7345582bec1SHong Zhang PetscFunctionReturn(0); 7355582bec1SHong Zhang } 736e14861a4SHong Zhang 737e14861a4SHong Zhang #undef __FUNCT__ 738eef31507SHong Zhang #define __FUNCT__ "MatWrapML_MPIAIJ" 739eef31507SHong Zhang PetscErrorCode MatWrapML_MPIAIJ(ML_Operator *mlmat,Mat *newmat) 740e14861a4SHong Zhang { 741ab718edeSHong Zhang struct ML_CSR_MSRdata *matdata = (struct ML_CSR_MSRdata *)mlmat->data; 742ab718edeSHong Zhang PetscInt *ml_cols=matdata->columns,*aj; 743ab718edeSHong Zhang PetscScalar *ml_vals=matdata->values,*aa; 744e14861a4SHong Zhang PetscErrorCode ierr; 745ab718edeSHong Zhang PetscInt i,j,k,*gordering; 7463e826d7bSSatish Balay PetscInt m=mlmat->outvec_leng,n,*nnzA,*nnzB,*nnz,nz_max,row; 747ab718edeSHong Zhang Mat A; 748e14861a4SHong Zhang 749e14861a4SHong Zhang PetscFunctionBegin; 750ab718edeSHong Zhang if (mlmat->getrow == NULL) SETERRQ(PETSC_ERR_ARG_NULL,"mlmat->getrow = NULL"); 751ab718edeSHong Zhang n = mlmat->invec_leng; 752e14861a4SHong Zhang if (m != n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"m %d must equal to n %d",m,n); 753e14861a4SHong Zhang 754f69a0ea3SMatthew Knepley ierr = MatCreate(mlmat->comm->USR_comm,&A);CHKERRQ(ierr); 755f69a0ea3SMatthew Knepley ierr = MatSetSizes(A,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 756ab718edeSHong Zhang ierr = MatSetType(A,MATMPIAIJ);CHKERRQ(ierr); 7573e826d7bSSatish Balay ierr = PetscMalloc3(m,PetscInt,&nnzA,m,PetscInt,&nnzB,m,PetscInt,&nnz);CHKERRQ(ierr); 7583e826d7bSSatish Balay 759e14861a4SHong Zhang nz_max = 0; 760e14861a4SHong Zhang for (i=0; i<m; i++){ 761ab718edeSHong Zhang nnz[i] = ml_cols[i+1] - ml_cols[i] + 1; 762e14861a4SHong Zhang if (nz_max < nnz[i]) nz_max = nnz[i]; 763ab718edeSHong Zhang nnzA[i] = 1; /* diag */ 764ab718edeSHong Zhang for (j=ml_cols[i]; j<ml_cols[i+1]; j++){ 765ab718edeSHong Zhang if (ml_cols[j] < m) nnzA[i]++; 766e14861a4SHong Zhang } 767e14861a4SHong Zhang nnzB[i] = nnz[i] - nnzA[i]; 768e14861a4SHong Zhang } 769ab718edeSHong Zhang ierr = MatMPIAIJSetPreallocation(A,0,nnzA,0,nnzB);CHKERRQ(ierr); 770e14861a4SHong Zhang 771ab718edeSHong Zhang /* insert mat values -- remap row and column indices */ 772ab718edeSHong Zhang nz_max++; 7737c307921SBarry Smith ierr = PetscMalloc(nz_max*(sizeof(PetscInt)+sizeof(PetscScalar)),&aj);CHKERRQ(ierr); 774ab718edeSHong Zhang aa = (PetscScalar*)(aj + nz_max); 7751e5ab15bSHong Zhang ML_build_global_numbering(mlmat,&gordering); 776e14861a4SHong Zhang for (i=0; i<m; i++){ 777e14861a4SHong Zhang row = gordering[i]; 778ab718edeSHong Zhang k = 0; 779ab718edeSHong Zhang /* diagonal entry */ 780ab718edeSHong Zhang aj[k] = row; aa[k++] = ml_vals[i]; 781ab718edeSHong Zhang /* off diagonal entries */ 782ab718edeSHong Zhang for (j=ml_cols[i]; j<ml_cols[i+1]; j++){ 783ab718edeSHong Zhang aj[k] = gordering[ml_cols[j]]; aa[k++] = ml_vals[j]; 784e14861a4SHong Zhang } 785ab718edeSHong Zhang ierr = MatSetValues(A,1,&row,nnz[i],aj,aa,INSERT_VALUES);CHKERRQ(ierr); 786e14861a4SHong Zhang } 787ab718edeSHong Zhang ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 788ab718edeSHong Zhang ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 789ab718edeSHong Zhang *newmat = A; 790e14861a4SHong Zhang 7913e826d7bSSatish Balay ierr = PetscFree3(nnzA,nnzB,nnz); 792ab718edeSHong Zhang ierr = PetscFree(aj);CHKERRQ(ierr); 793e14861a4SHong Zhang PetscFunctionReturn(0); 794e14861a4SHong Zhang } 795