1dba47a55SKris Buschelman #define PETSCKSP_DLL 2ab718edeSHong Zhang 35582bec1SHong Zhang /* 4759e7b9cSHong Zhang Provides an interface to the ML 5.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 #include <math.h> 122cf39c26SSatish Balay EXTERN_C_BEGIN 132cf39c26SSatish Balay #include "ml_config.h" 145582bec1SHong Zhang #include "ml_include.h" 155582bec1SHong Zhang EXTERN_C_END 165582bec1SHong Zhang 175582bec1SHong Zhang /* The context (data structure) at each grid level */ 185582bec1SHong Zhang typedef struct { 195582bec1SHong Zhang Vec x,b,r; /* global vectors */ 205582bec1SHong Zhang Mat A,P,R; 215582bec1SHong Zhang KSP ksp; 225582bec1SHong Zhang } GridCtx; 235582bec1SHong Zhang 245582bec1SHong Zhang /* The context used to input PETSc matrix into ML at fine grid */ 255582bec1SHong Zhang typedef struct { 26573998d7SHong Zhang Mat A; /* Petsc matrix in aij format */ 27573998d7SHong Zhang Mat Aloc; /* local portion of A to be used by ML */ 2824a42b14SHong Zhang Vec x,y; 295582bec1SHong Zhang ML_Operator *mlmat; 305582bec1SHong Zhang PetscScalar *pwork; /* tmp array used by PetscML_comm() */ 315582bec1SHong Zhang } FineGridCtx; 325582bec1SHong Zhang 335582bec1SHong Zhang /* The context associates a ML matrix with a PETSc shell matrix */ 345582bec1SHong Zhang typedef struct { 355582bec1SHong Zhang Mat A; /* PETSc shell matrix associated with mlmat */ 365582bec1SHong Zhang ML_Operator *mlmat; /* ML matrix assorciated with A */ 375582bec1SHong Zhang Vec y; 385582bec1SHong Zhang } Mat_MLShell; 395582bec1SHong Zhang 405582bec1SHong Zhang /* Private context for the ML preconditioner */ 415582bec1SHong Zhang typedef struct { 425582bec1SHong Zhang ML *ml_object; 435582bec1SHong Zhang ML_Aggregate *agg_object; 445582bec1SHong Zhang GridCtx *gridctx; 455582bec1SHong Zhang FineGridCtx *PetscMLdata; 46573998d7SHong Zhang PetscInt Nlevels,MaxNlevels,MaxCoarseSize,CoarsenScheme; 475582bec1SHong Zhang PetscReal Threshold,DampingFactor; 485582bec1SHong Zhang PetscTruth SpectralNormScheme_Anorm; 49573998d7SHong Zhang PetscMPIInt size; /* size of communicator for pc->pmat */ 505582bec1SHong Zhang PetscErrorCode (*PCSetUp)(PC); 515582bec1SHong Zhang PetscErrorCode (*PCDestroy)(PC); 525582bec1SHong Zhang } PC_ML; 5341ca0015SHong Zhang 5441ca0015SHong Zhang extern int PetscML_getrow(ML_Operator *ML_data,int N_requested_rows,int requested_rows[], 555582bec1SHong Zhang int allocated_space,int columns[],double values[],int row_lengths[]); 5641ca0015SHong Zhang extern int PetscML_matvec(ML_Operator *ML_data, int in_length, double p[], int out_length,double ap[]); 575582bec1SHong Zhang extern int PetscML_comm(double x[], void *ML_data); 585582bec1SHong Zhang extern PetscErrorCode MatMult_ML(Mat,Vec,Vec); 595582bec1SHong Zhang extern PetscErrorCode MatMultAdd_ML(Mat,Vec,Vec,Vec); 6075179d2cSHong Zhang extern PetscErrorCode MatConvert_MPIAIJ_ML(Mat,MatType,MatReuse,Mat*); 615582bec1SHong Zhang extern PetscErrorCode MatDestroy_ML(Mat); 62573998d7SHong Zhang extern PetscErrorCode MatWrapML_SeqAIJ(ML_Operator*,MatReuse,Mat*); 63eef31507SHong Zhang extern PetscErrorCode MatWrapML_MPIAIJ(ML_Operator*,Mat*); 64573998d7SHong Zhang extern PetscErrorCode MatWrapML_SHELL(ML_Operator*,MatReuse,Mat*); 65573998d7SHong Zhang extern PetscErrorCode PetscContainerDestroy_PC_ML(void *); 665582bec1SHong Zhang 675582bec1SHong Zhang /* -------------------------------------------------------------------------- */ 685582bec1SHong Zhang /* 695582bec1SHong Zhang PCSetUp_ML - Prepares for the use of the ML preconditioner 705582bec1SHong Zhang by setting data structures and options. 715582bec1SHong Zhang 725582bec1SHong Zhang Input Parameter: 735582bec1SHong Zhang . pc - the preconditioner context 745582bec1SHong Zhang 755582bec1SHong Zhang Application Interface Routine: PCSetUp() 765582bec1SHong Zhang 775582bec1SHong Zhang Notes: 785582bec1SHong Zhang The interface routine PCSetUp() is not usually called directly by 795582bec1SHong Zhang the user, but instead is called by PCApply() if necessary. 805582bec1SHong Zhang */ 816ca4d86aSHong Zhang extern PetscErrorCode PCSetFromOptions_MG(PC); 825582bec1SHong Zhang #undef __FUNCT__ 835582bec1SHong Zhang #define __FUNCT__ "PCSetUp_ML" 846ca4d86aSHong Zhang PetscErrorCode PCSetUp_ML(PC pc) 855582bec1SHong Zhang { 865582bec1SHong Zhang PetscErrorCode ierr; 87eef31507SHong Zhang PetscMPIInt size; 885582bec1SHong Zhang FineGridCtx *PetscMLdata; 895582bec1SHong Zhang ML *ml_object; 905582bec1SHong Zhang ML_Aggregate *agg_object; 915582bec1SHong Zhang ML_Operator *mlmat; 92ac346b81SHong Zhang PetscInt nlocal_allcols,Nlevels,mllevel,level,level1,m,fine_level; 935582bec1SHong Zhang Mat A,Aloc; 945582bec1SHong Zhang GridCtx *gridctx; 955582bec1SHong Zhang PC_ML *pc_ml=PETSC_NULL; 96776b82aeSLisandro Dalcin PetscContainer container; 97573998d7SHong Zhang MatReuse reuse = MAT_INITIAL_MATRIX; 985582bec1SHong Zhang 995582bec1SHong Zhang PetscFunctionBegin; 1005582bec1SHong Zhang ierr = PetscObjectQuery((PetscObject)pc,"PC_ML",(PetscObject *)&container);CHKERRQ(ierr); 1015582bec1SHong Zhang if (container) { 102776b82aeSLisandro Dalcin ierr = PetscContainerGetPointer(container,(void **)&pc_ml);CHKERRQ(ierr); 1035582bec1SHong Zhang } else { 1045582bec1SHong Zhang SETERRQ(PETSC_ERR_ARG_NULL,"Container does not exit"); 1055582bec1SHong Zhang } 1065582bec1SHong Zhang 107573998d7SHong Zhang if (pc->setupcalled){ 108573998d7SHong Zhang if (pc->flag == SAME_NONZERO_PATTERN){ 109573998d7SHong Zhang reuse = MAT_REUSE_MATRIX; 110573998d7SHong Zhang PetscMLdata = pc_ml->PetscMLdata; 111573998d7SHong Zhang gridctx = pc_ml->gridctx; 112573998d7SHong Zhang /* ML objects cannot be reused */ 113573998d7SHong Zhang ML_Destroy(&pc_ml->ml_object); 114573998d7SHong Zhang ML_Aggregate_Destroy(&pc_ml->agg_object); 115573998d7SHong Zhang } else { 116573998d7SHong Zhang PC_ML *pc_ml_new = PETSC_NULL; 117573998d7SHong Zhang PetscContainer container_new; 11838f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_ML,&pc_ml_new);CHKERRQ(ierr); 119573998d7SHong Zhang ierr = PetscContainerCreate(PETSC_COMM_SELF,&container_new);CHKERRQ(ierr); 120573998d7SHong Zhang ierr = PetscContainerSetPointer(container_new,pc_ml_new);CHKERRQ(ierr); 121573998d7SHong Zhang ierr = PetscContainerSetUserDestroy(container_new,PetscContainerDestroy_PC_ML);CHKERRQ(ierr); 122573998d7SHong Zhang ierr = PetscObjectCompose((PetscObject)pc,"PC_ML",(PetscObject)container_new);CHKERRQ(ierr); 123573998d7SHong Zhang 124573998d7SHong Zhang ierr = PetscMemcpy(pc_ml_new,pc_ml,sizeof(PC_ML));CHKERRQ(ierr); 125573998d7SHong Zhang ierr = PetscContainerDestroy(container);CHKERRQ(ierr); 126573998d7SHong Zhang pc_ml = pc_ml_new; 127573998d7SHong Zhang } 128573998d7SHong Zhang } 129573998d7SHong Zhang 1305582bec1SHong Zhang /* setup special features of PCML */ 1315582bec1SHong Zhang /*--------------------------------*/ 1325582bec1SHong Zhang /* covert A to Aloc to be used by ML at fine grid */ 1335582bec1SHong Zhang A = pc->pmat; 1347adad957SLisandro Dalcin ierr = MPI_Comm_size(((PetscObject)A)->comm,&size);CHKERRQ(ierr); 1355582bec1SHong Zhang pc_ml->size = size; 1365582bec1SHong Zhang if (size > 1){ 137573998d7SHong Zhang if (reuse) Aloc = PetscMLdata->Aloc; 138573998d7SHong Zhang ierr = MatConvert_MPIAIJ_ML(A,PETSC_NULL,reuse,&Aloc);CHKERRQ(ierr); 1395582bec1SHong Zhang } else { 1405582bec1SHong Zhang Aloc = A; 1415582bec1SHong Zhang } 1425582bec1SHong Zhang 1435582bec1SHong Zhang /* create and initialize struct 'PetscMLdata' */ 144573998d7SHong Zhang if (!reuse){ 14538f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,FineGridCtx,&PetscMLdata);CHKERRQ(ierr); 1465582bec1SHong Zhang pc_ml->PetscMLdata = PetscMLdata; 147573998d7SHong Zhang ierr = PetscMalloc((Aloc->cmap.n+1)*sizeof(PetscScalar),&PetscMLdata->pwork);CHKERRQ(ierr); 1485582bec1SHong Zhang 14924a42b14SHong Zhang ierr = VecCreate(PETSC_COMM_SELF,&PetscMLdata->x);CHKERRQ(ierr); 150a803a431SHong Zhang ierr = VecSetSizes(PetscMLdata->x,Aloc->cmap.n,Aloc->cmap.n);CHKERRQ(ierr); 15124a42b14SHong Zhang ierr = VecSetType(PetscMLdata->x,VECSEQ);CHKERRQ(ierr); 15224a42b14SHong Zhang 15324a42b14SHong Zhang ierr = VecCreate(PETSC_COMM_SELF,&PetscMLdata->y);CHKERRQ(ierr); 154a803a431SHong Zhang ierr = VecSetSizes(PetscMLdata->y,A->rmap.n,PETSC_DECIDE);CHKERRQ(ierr); 15524a42b14SHong Zhang ierr = VecSetType(PetscMLdata->y,VECSEQ);CHKERRQ(ierr); 156573998d7SHong Zhang } 157573998d7SHong Zhang PetscMLdata->A = A; 158573998d7SHong Zhang PetscMLdata->Aloc = Aloc; 15924a42b14SHong Zhang 1605582bec1SHong Zhang /* create ML discretization matrix at fine grid */ 16145cf47abSHong Zhang /* ML requires input of fine-grid matrix. It determines nlevels. */ 1625582bec1SHong Zhang ierr = MatGetSize(Aloc,&m,&nlocal_allcols);CHKERRQ(ierr); 1635582bec1SHong Zhang ML_Create(&ml_object,pc_ml->MaxNlevels); 164573998d7SHong Zhang pc_ml->ml_object = ml_object; 1655582bec1SHong Zhang ML_Init_Amatrix(ml_object,0,m,m,PetscMLdata); 1665582bec1SHong Zhang ML_Set_Amatrix_Getrow(ml_object,0,PetscML_getrow,PetscML_comm,nlocal_allcols); 1675582bec1SHong Zhang ML_Set_Amatrix_Matvec(ml_object,0,PetscML_matvec); 1685582bec1SHong Zhang 1695582bec1SHong Zhang /* aggregation */ 1705582bec1SHong Zhang ML_Aggregate_Create(&agg_object); 171573998d7SHong Zhang pc_ml->agg_object = agg_object; 172573998d7SHong Zhang 1735582bec1SHong Zhang ML_Aggregate_Set_MaxCoarseSize(agg_object,pc_ml->MaxCoarseSize); 1745582bec1SHong Zhang /* set options */ 1755582bec1SHong Zhang switch (pc_ml->CoarsenScheme) { 1765582bec1SHong Zhang case 1: 1775582bec1SHong Zhang ML_Aggregate_Set_CoarsenScheme_Coupled(agg_object);break; 1785582bec1SHong Zhang case 2: 1795582bec1SHong Zhang ML_Aggregate_Set_CoarsenScheme_MIS(agg_object);break; 1805582bec1SHong Zhang case 3: 1815582bec1SHong Zhang ML_Aggregate_Set_CoarsenScheme_METIS(agg_object);break; 1825582bec1SHong Zhang } 1835582bec1SHong Zhang ML_Aggregate_Set_Threshold(agg_object,pc_ml->Threshold); 1845582bec1SHong Zhang ML_Aggregate_Set_DampingFactor(agg_object,pc_ml->DampingFactor); 1855582bec1SHong Zhang if (pc_ml->SpectralNormScheme_Anorm){ 1865582bec1SHong Zhang ML_Aggregate_Set_SpectralNormScheme_Anorm(agg_object); 1875582bec1SHong Zhang } 1885582bec1SHong Zhang 1895582bec1SHong Zhang Nlevels = ML_Gen_MGHierarchy_UsingAggregation(ml_object,0,ML_INCREASING,agg_object); 1905582bec1SHong Zhang if (Nlevels<=0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Nlevels %d must > 0",Nlevels); 191573998d7SHong Zhang if (pc->setupcalled && pc_ml->Nlevels != Nlevels) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"previous Nlevels %D and current Nlevels %d must be same", pc_ml->Nlevels,Nlevels); 192573998d7SHong Zhang pc_ml->Nlevels = Nlevels; 193*aa85bbbfSHong Zhang fine_level = Nlevels - 1; 194573998d7SHong Zhang if (!pc->setupcalled){ 19597177400SBarry Smith ierr = PCMGSetLevels(pc,Nlevels,PETSC_NULL);CHKERRQ(ierr); 196*aa85bbbfSHong Zhang /* set default smoothers */ 197*aa85bbbfSHong Zhang KSP smoother; 198*aa85bbbfSHong Zhang PC subpc; 199*aa85bbbfSHong Zhang for (level=1; level<=fine_level; level++){ 200*aa85bbbfSHong Zhang if (size == 1){ 201*aa85bbbfSHong Zhang ierr = PCMGGetSmoother(pc,level,&smoother);CHKERRQ(ierr); 202*aa85bbbfSHong Zhang ierr = KSPSetType(smoother,KSPRICHARDSON);CHKERRQ(ierr); 203*aa85bbbfSHong Zhang ierr = KSPGetPC(smoother,&subpc);CHKERRQ(ierr); 204*aa85bbbfSHong Zhang ierr = PCSetType(subpc,PCSOR);CHKERRQ(ierr); 205*aa85bbbfSHong Zhang } else { 206*aa85bbbfSHong Zhang ierr = PCMGGetSmoother(pc,level,&smoother);CHKERRQ(ierr); 207*aa85bbbfSHong Zhang ierr = KSPSetType(smoother,KSPRICHARDSON);CHKERRQ(ierr); 208*aa85bbbfSHong Zhang ierr = KSPGetPC(smoother,&subpc);CHKERRQ(ierr); 209*aa85bbbfSHong Zhang ierr = PCSetType(subpc,PCSOR);CHKERRQ(ierr); 210*aa85bbbfSHong Zhang } 211*aa85bbbfSHong Zhang } 21297177400SBarry Smith ierr = PCSetFromOptions_MG(pc);CHKERRQ(ierr); /* should be called in PCSetFromOptions_ML(), but cannot be called prior to PCMGSetLevels() */ 213573998d7SHong Zhang } 2145582bec1SHong Zhang 215573998d7SHong Zhang if (!reuse){ 2165582bec1SHong Zhang ierr = PetscMalloc(Nlevels*sizeof(GridCtx),&gridctx);CHKERRQ(ierr); 2175582bec1SHong Zhang pc_ml->gridctx = gridctx; 218573998d7SHong Zhang } 2195582bec1SHong Zhang 2205582bec1SHong Zhang /* wrap ML matrices by PETSc shell matrices at coarsened grids. 2215582bec1SHong Zhang Level 0 is the finest grid for ML, but coarsest for PETSc! */ 222e14861a4SHong Zhang gridctx[fine_level].A = A; 223573998d7SHong Zhang 224e14861a4SHong Zhang level = fine_level - 1; 225ab718edeSHong Zhang if (size == 1){ /* convert ML P, R and A into seqaij format */ 2265582bec1SHong Zhang for (mllevel=1; mllevel<Nlevels; mllevel++){ 227e14861a4SHong Zhang mlmat = &(ml_object->Pmat[mllevel]); 228573998d7SHong Zhang ierr = MatWrapML_SeqAIJ(mlmat,reuse,&gridctx[level].P);CHKERRQ(ierr); 229e14861a4SHong Zhang mlmat = &(ml_object->Rmat[mllevel-1]); 230573998d7SHong Zhang ierr = MatWrapML_SeqAIJ(mlmat,reuse,&gridctx[level].R);CHKERRQ(ierr); 231573998d7SHong Zhang 232573998d7SHong Zhang mlmat = &(ml_object->Amat[mllevel]); 233573998d7SHong Zhang if (reuse){ 234573998d7SHong Zhang /* ML matrix A changes sparse pattern although PETSc A doesn't, thus gridctx[level].A must be recreated! */ 235573998d7SHong Zhang ierr = MatDestroy(gridctx[level].A);CHKERRQ(ierr); 236573998d7SHong Zhang } 237573998d7SHong Zhang ierr = MatWrapML_SeqAIJ(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].A);CHKERRQ(ierr); 2385582bec1SHong Zhang level--; 2395582bec1SHong Zhang } 240ab718edeSHong Zhang } else { /* convert ML P and R into shell format, ML A into mpiaij format */ 2415582bec1SHong Zhang for (mllevel=1; mllevel<Nlevels; mllevel++){ 2425582bec1SHong Zhang mlmat = &(ml_object->Pmat[mllevel]); 243573998d7SHong Zhang ierr = MatWrapML_SHELL(mlmat,reuse,&gridctx[level].P);CHKERRQ(ierr); 244ab718edeSHong Zhang mlmat = &(ml_object->Rmat[mllevel-1]); 245573998d7SHong Zhang ierr = MatWrapML_SHELL(mlmat,reuse,&gridctx[level].R);CHKERRQ(ierr); 246573998d7SHong Zhang 2475582bec1SHong Zhang mlmat = &(ml_object->Amat[mllevel]); 248573998d7SHong Zhang if (reuse){ 249573998d7SHong Zhang ierr = MatDestroy(gridctx[level].A);CHKERRQ(ierr); 250573998d7SHong Zhang } 251eef31507SHong Zhang ierr = MatWrapML_MPIAIJ(mlmat,&gridctx[level].A);CHKERRQ(ierr); 2525582bec1SHong Zhang level--; 2535582bec1SHong Zhang } 2545582bec1SHong Zhang } 2555582bec1SHong Zhang 256573998d7SHong Zhang /* create vectors and ksp at all levels */ 257573998d7SHong Zhang if (!reuse){ 258ac346b81SHong Zhang for (level=0; level<fine_level; level++){ 259573998d7SHong Zhang level1 = level + 1; 260e64afeacSLisandro Dalcin ierr = VecCreate(((PetscObject)gridctx[level].A)->comm,&gridctx[level].x);CHKERRQ(ierr); 261a803a431SHong Zhang ierr = VecSetSizes(gridctx[level].x,gridctx[level].A->cmap.n,PETSC_DECIDE);CHKERRQ(ierr); 2625582bec1SHong Zhang ierr = VecSetType(gridctx[level].x,VECMPI);CHKERRQ(ierr); 26397177400SBarry Smith ierr = PCMGSetX(pc,level,gridctx[level].x);CHKERRQ(ierr); 2645582bec1SHong Zhang 265e64afeacSLisandro Dalcin ierr = VecCreate(((PetscObject)gridctx[level].A)->comm,&gridctx[level].b);CHKERRQ(ierr); 266a803a431SHong Zhang ierr = VecSetSizes(gridctx[level].b,gridctx[level].A->rmap.n,PETSC_DECIDE);CHKERRQ(ierr); 2675582bec1SHong Zhang ierr = VecSetType(gridctx[level].b,VECMPI);CHKERRQ(ierr); 26897177400SBarry Smith ierr = PCMGSetRhs(pc,level,gridctx[level].b);CHKERRQ(ierr); 269ac346b81SHong Zhang 270e64afeacSLisandro Dalcin ierr = VecCreate(((PetscObject)gridctx[level1].A)->comm,&gridctx[level1].r);CHKERRQ(ierr); 271a803a431SHong Zhang ierr = VecSetSizes(gridctx[level1].r,gridctx[level1].A->rmap.n,PETSC_DECIDE);CHKERRQ(ierr); 272ac346b81SHong Zhang ierr = VecSetType(gridctx[level1].r,VECMPI);CHKERRQ(ierr); 27397177400SBarry Smith ierr = PCMGSetR(pc,level1,gridctx[level1].r);CHKERRQ(ierr); 274ac346b81SHong Zhang 2755582bec1SHong Zhang if (level == 0){ 27697177400SBarry Smith ierr = PCMGGetCoarseSolve(pc,&gridctx[level].ksp);CHKERRQ(ierr); 2775582bec1SHong Zhang } else { 27897177400SBarry Smith ierr = PCMGGetSmoother(pc,level,&gridctx[level].ksp);CHKERRQ(ierr); 279573998d7SHong Zhang } 280573998d7SHong Zhang } 281573998d7SHong Zhang ierr = PCMGGetSmoother(pc,fine_level,&gridctx[fine_level].ksp);CHKERRQ(ierr); 282573998d7SHong Zhang } 283573998d7SHong Zhang 284573998d7SHong Zhang /* create coarse level and the interpolation between the levels */ 285573998d7SHong Zhang for (level=0; level<fine_level; level++){ 286573998d7SHong Zhang level1 = level + 1; 287aea2a34eSBarry Smith ierr = PCMGSetInterpolation(pc,level1,gridctx[level].P);CHKERRQ(ierr); 288573998d7SHong Zhang ierr = PCMGSetRestriction(pc,level1,gridctx[level].R);CHKERRQ(ierr); 289573998d7SHong Zhang if (level > 0){ 29097177400SBarry Smith ierr = PCMGSetResidual(pc,level,PCMGDefaultResidual,gridctx[level].A);CHKERRQ(ierr); 2915582bec1SHong Zhang } 2925582bec1SHong Zhang ierr = KSPSetOperators(gridctx[level].ksp,gridctx[level].A,gridctx[level].A,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 2935582bec1SHong Zhang } 29497177400SBarry Smith ierr = PCMGSetResidual(pc,fine_level,PCMGDefaultResidual,gridctx[fine_level].A);CHKERRQ(ierr); 295ac346b81SHong Zhang ierr = KSPSetOperators(gridctx[fine_level].ksp,gridctx[level].A,gridctx[fine_level].A,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 2965582bec1SHong Zhang 2975582bec1SHong Zhang /* now call PCSetUp_MG() */ 298573998d7SHong Zhang /*-------------------------------*/ 2995582bec1SHong Zhang ierr = (*pc_ml->PCSetUp)(pc);CHKERRQ(ierr); 3005582bec1SHong Zhang PetscFunctionReturn(0); 3015582bec1SHong Zhang } 3025582bec1SHong Zhang 3035582bec1SHong Zhang #undef __FUNCT__ 304776b82aeSLisandro Dalcin #define __FUNCT__ "PetscContainerDestroy_PC_ML" 305776b82aeSLisandro Dalcin PetscErrorCode PetscContainerDestroy_PC_ML(void *ptr) 3065582bec1SHong Zhang { 3075582bec1SHong Zhang PetscErrorCode ierr; 3085582bec1SHong Zhang PC_ML *pc_ml = (PC_ML*)ptr; 309573998d7SHong Zhang PetscInt level,fine_level=pc_ml->Nlevels-1; 3105582bec1SHong Zhang 3115582bec1SHong Zhang PetscFunctionBegin; 3125582bec1SHong Zhang ML_Aggregate_Destroy(&pc_ml->agg_object); 3135582bec1SHong Zhang ML_Destroy(&pc_ml->ml_object); 3145582bec1SHong Zhang 31538f2d2fdSLisandro Dalcin if (pc_ml->PetscMLdata) { 3165582bec1SHong Zhang ierr = PetscFree(pc_ml->PetscMLdata->pwork);CHKERRQ(ierr); 31738f2d2fdSLisandro Dalcin if (pc_ml->size > 1) {ierr = MatDestroy(pc_ml->PetscMLdata->Aloc);CHKERRQ(ierr);} 318ac346b81SHong Zhang if (pc_ml->PetscMLdata->x){ierr = VecDestroy(pc_ml->PetscMLdata->x);CHKERRQ(ierr);} 319ac346b81SHong Zhang if (pc_ml->PetscMLdata->y){ierr = VecDestroy(pc_ml->PetscMLdata->y);CHKERRQ(ierr);} 32038f2d2fdSLisandro Dalcin } 3215582bec1SHong Zhang ierr = PetscFree(pc_ml->PetscMLdata);CHKERRQ(ierr); 3225582bec1SHong Zhang 323573998d7SHong Zhang for (level=0; level<fine_level; level++){ 324ac346b81SHong Zhang if (pc_ml->gridctx[level].A){ierr = MatDestroy(pc_ml->gridctx[level].A);CHKERRQ(ierr);} 325ac346b81SHong Zhang if (pc_ml->gridctx[level].P){ierr = MatDestroy(pc_ml->gridctx[level].P);CHKERRQ(ierr);} 326ac346b81SHong Zhang if (pc_ml->gridctx[level].R){ierr = MatDestroy(pc_ml->gridctx[level].R);CHKERRQ(ierr);} 327ac346b81SHong Zhang if (pc_ml->gridctx[level].x){ierr = VecDestroy(pc_ml->gridctx[level].x);CHKERRQ(ierr);} 328ac346b81SHong Zhang if (pc_ml->gridctx[level].b){ierr = VecDestroy(pc_ml->gridctx[level].b);CHKERRQ(ierr);} 329ac346b81SHong Zhang if (pc_ml->gridctx[level+1].r){ierr = VecDestroy(pc_ml->gridctx[level+1].r);CHKERRQ(ierr);} 3305582bec1SHong Zhang } 3315582bec1SHong Zhang ierr = PetscFree(pc_ml->gridctx);CHKERRQ(ierr); 3325582bec1SHong Zhang ierr = PetscFree(pc_ml);CHKERRQ(ierr); 3335582bec1SHong Zhang PetscFunctionReturn(0); 3345582bec1SHong Zhang } 3355582bec1SHong Zhang /* -------------------------------------------------------------------------- */ 3365582bec1SHong Zhang /* 3375582bec1SHong Zhang PCDestroy_ML - Destroys the private context for the ML preconditioner 3385582bec1SHong Zhang that was created with PCCreate_ML(). 3395582bec1SHong Zhang 3405582bec1SHong Zhang Input Parameter: 3415582bec1SHong Zhang . pc - the preconditioner context 3425582bec1SHong Zhang 3435582bec1SHong Zhang Application Interface Routine: PCDestroy() 3445582bec1SHong Zhang */ 3455582bec1SHong Zhang #undef __FUNCT__ 3465582bec1SHong Zhang #define __FUNCT__ "PCDestroy_ML" 3476ca4d86aSHong Zhang PetscErrorCode PCDestroy_ML(PC pc) 3485582bec1SHong Zhang { 3495582bec1SHong Zhang PetscErrorCode ierr; 3505582bec1SHong Zhang PC_ML *pc_ml=PETSC_NULL; 351776b82aeSLisandro Dalcin PetscContainer container; 3525582bec1SHong Zhang 3535582bec1SHong Zhang PetscFunctionBegin; 3545582bec1SHong Zhang ierr = PetscObjectQuery((PetscObject)pc,"PC_ML",(PetscObject *)&container);CHKERRQ(ierr); 3555582bec1SHong Zhang if (container) { 356776b82aeSLisandro Dalcin ierr = PetscContainerGetPointer(container,(void **)&pc_ml);CHKERRQ(ierr); 3575582bec1SHong Zhang pc->ops->destroy = pc_ml->PCDestroy; 3585582bec1SHong Zhang } else { 3595582bec1SHong Zhang SETERRQ(PETSC_ERR_ARG_NULL,"Container does not exit"); 3605582bec1SHong Zhang } 3619cb0ec93SHong Zhang /* detach pc and PC_ML and dereference container */ 36238f2d2fdSLisandro Dalcin ierr = PetscContainerDestroy(container);CHKERRQ(ierr); 3635582bec1SHong Zhang ierr = PetscObjectCompose((PetscObject)pc,"PC_ML",0);CHKERRQ(ierr); 36438f2d2fdSLisandro Dalcin if (pc->ops->destroy) { 3655582bec1SHong Zhang ierr = (*pc->ops->destroy)(pc);CHKERRQ(ierr); 36638f2d2fdSLisandro Dalcin } 3675582bec1SHong Zhang PetscFunctionReturn(0); 3685582bec1SHong Zhang } 3695582bec1SHong Zhang 3705582bec1SHong Zhang #undef __FUNCT__ 3715582bec1SHong Zhang #define __FUNCT__ "PCSetFromOptions_ML" 3726ca4d86aSHong Zhang PetscErrorCode PCSetFromOptions_ML(PC pc) 3735582bec1SHong Zhang { 3745582bec1SHong Zhang PetscErrorCode ierr; 37538f2d2fdSLisandro Dalcin PetscInt indx,m,PrintLevel; 3765582bec1SHong Zhang PetscTruth flg; 3775582bec1SHong Zhang const char *scheme[] = {"Uncoupled","Coupled","MIS","METIS"}; 3785582bec1SHong Zhang PC_ML *pc_ml=PETSC_NULL; 379776b82aeSLisandro Dalcin PetscContainer container; 3809dcbbd2bSBarry Smith PCMGType mgtype; 3815582bec1SHong Zhang 3825582bec1SHong Zhang PetscFunctionBegin; 3835582bec1SHong Zhang ierr = PetscObjectQuery((PetscObject)pc,"PC_ML",(PetscObject *)&container);CHKERRQ(ierr); 3845582bec1SHong Zhang if (container) { 385776b82aeSLisandro Dalcin ierr = PetscContainerGetPointer(container,(void **)&pc_ml);CHKERRQ(ierr); 3865582bec1SHong Zhang } else { 3875582bec1SHong Zhang SETERRQ(PETSC_ERR_ARG_NULL,"Container does not exit"); 3885582bec1SHong Zhang } 3896ca4d86aSHong Zhang 3905582bec1SHong Zhang /* inherited MG options */ 3916ca4d86aSHong Zhang ierr = PetscOptionsHead("Multigrid options(inherited)");CHKERRQ(ierr); 3925582bec1SHong Zhang ierr = PetscOptionsInt("-pc_mg_cycles","1 for V cycle, 2 for W-cycle","MGSetCycles",1,&m,&flg);CHKERRQ(ierr); 3935582bec1SHong Zhang ierr = PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","MGSetNumberSmoothUp",1,&m,&flg);CHKERRQ(ierr); 3945582bec1SHong Zhang ierr = PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","MGSetNumberSmoothDown",1,&m,&flg);CHKERRQ(ierr); 3959dcbbd2bSBarry Smith ierr = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)PC_MG_MULTIPLICATIVE,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr); 3965582bec1SHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 3975582bec1SHong Zhang 3985582bec1SHong Zhang /* ML options */ 3995582bec1SHong Zhang ierr = PetscOptionsHead("ML options");CHKERRQ(ierr); 4005582bec1SHong Zhang /* set defaults */ 4015582bec1SHong Zhang PrintLevel = 0; 4025582bec1SHong Zhang indx = 0; 4035582bec1SHong Zhang ierr = PetscOptionsInt("-pc_ml_PrintLevel","Print level","ML_Set_PrintLevel",PrintLevel,&PrintLevel,PETSC_NULL);CHKERRQ(ierr); 4045582bec1SHong Zhang ML_Set_PrintLevel(PrintLevel); 405573998d7SHong Zhang ierr = PetscOptionsInt("-pc_ml_maxNlevels","Maximum number of levels","None",pc_ml->MaxNlevels,&pc_ml->MaxNlevels,PETSC_NULL);CHKERRQ(ierr); 406573998d7SHong Zhang ierr = PetscOptionsInt("-pc_ml_maxCoarseSize","Maximum coarsest mesh size","ML_Aggregate_Set_MaxCoarseSize",pc_ml->MaxCoarseSize,&pc_ml->MaxCoarseSize,PETSC_NULL);CHKERRQ(ierr); 407573998d7SHong Zhang ierr = PetscOptionsEList("-pc_ml_CoarsenScheme","Aggregate Coarsen Scheme","ML_Aggregate_Set_CoarsenScheme_*",scheme,4,scheme[0],&indx,PETSC_NULL);CHKERRQ(ierr); /* ??? */ 4085582bec1SHong Zhang pc_ml->CoarsenScheme = indx; 4095582bec1SHong Zhang 410573998d7SHong Zhang ierr = PetscOptionsReal("-pc_ml_DampingFactor","P damping factor","ML_Aggregate_Set_DampingFactor",pc_ml->DampingFactor,&pc_ml->DampingFactor,PETSC_NULL);CHKERRQ(ierr); 4115582bec1SHong Zhang 412573998d7SHong Zhang ierr = PetscOptionsReal("-pc_ml_Threshold","Smoother drop tol","ML_Aggregate_Set_Threshold",pc_ml->Threshold,&pc_ml->Threshold,PETSC_NULL);CHKERRQ(ierr); 4135582bec1SHong Zhang 414573998d7SHong Zhang ierr = PetscOptionsTruth("-pc_ml_SpectralNormScheme_Anorm","Method used for estimating spectral radius","ML_Aggregate_Set_SpectralNormScheme_Anorm",pc_ml->SpectralNormScheme_Anorm,&pc_ml->SpectralNormScheme_Anorm,PETSC_NULL); 4155582bec1SHong Zhang 4165582bec1SHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 4175582bec1SHong Zhang PetscFunctionReturn(0); 4185582bec1SHong Zhang } 4195582bec1SHong Zhang 4205582bec1SHong Zhang /* -------------------------------------------------------------------------- */ 4215582bec1SHong Zhang /* 4225582bec1SHong Zhang PCCreate_ML - Creates a ML preconditioner context, PC_ML, 4235582bec1SHong Zhang and sets this as the private data within the generic preconditioning 4245582bec1SHong Zhang context, PC, that was created within PCCreate(). 4255582bec1SHong Zhang 4265582bec1SHong Zhang Input Parameter: 4275582bec1SHong Zhang . pc - the preconditioner context 4285582bec1SHong Zhang 4295582bec1SHong Zhang Application Interface Routine: PCCreate() 4305582bec1SHong Zhang */ 4315582bec1SHong Zhang 4325582bec1SHong Zhang /*MC 4331e5ab15bSHong Zhang PCML - Use algebraic multigrid preconditioning. This preconditioner requires you provide 4345582bec1SHong Zhang fine grid discretization matrix. The coarser grid matrices and restriction/interpolation 4356ca4d86aSHong Zhang operators are computed by ML, with the matrices coverted to PETSc matrices in aij format 4366ca4d86aSHong Zhang and the restriction/interpolation operators wrapped as PETSc shell matrices. 4375582bec1SHong Zhang 4386ca4d86aSHong Zhang Options Database Key: 4396ca4d86aSHong Zhang Multigrid options(inherited) 4406ca4d86aSHong Zhang + -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (MGSetCycles) 4416ca4d86aSHong Zhang . -pc_mg_smoothup <1>: Number of post-smoothing steps (MGSetNumberSmoothUp) 4426ca4d86aSHong Zhang . -pc_mg_smoothdown <1>: Number of pre-smoothing steps (MGSetNumberSmoothDown) 443f41ab451SVictor Eijkhout - -pc_mg_type <multiplicative>: (one of) additive multiplicative full cascade kascade 4446ca4d86aSHong Zhang 44551f519a2SBarry Smith ML options: 4466ca4d86aSHong Zhang + -pc_ml_PrintLevel <0>: Print level (ML_Set_PrintLevel) 4476ca4d86aSHong Zhang . -pc_ml_maxNlevels <10>: Maximum number of levels (None) 4486ca4d86aSHong Zhang . -pc_ml_maxCoarseSize <1>: Maximum coarsest mesh size (ML_Aggregate_Set_MaxCoarseSize) 449f41ab451SVictor Eijkhout . -pc_ml_CoarsenScheme <Uncoupled>: (one of) Uncoupled Coupled MIS METIS 4506ca4d86aSHong Zhang . -pc_ml_DampingFactor <1.33333>: P damping factor (ML_Aggregate_Set_DampingFactor) 4516ca4d86aSHong Zhang . -pc_ml_Threshold <0>: Smoother drop tol (ML_Aggregate_Set_Threshold) 452f41ab451SVictor Eijkhout - -pc_ml_SpectralNormScheme_Anorm <false>: Method used for estimating spectral radius (ML_Aggregate_Set_SpectralNormScheme_Anorm) 4535582bec1SHong Zhang 4545582bec1SHong Zhang Level: intermediate 4555582bec1SHong Zhang 4565582bec1SHong Zhang Concepts: multigrid 4575582bec1SHong Zhang 4585582bec1SHong Zhang .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, 45997177400SBarry Smith PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), MPSetCycles(), PCMGSetNumberSmoothDown(), 46097177400SBarry Smith PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(), 46197177400SBarry Smith PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(), 46297177400SBarry Smith PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR() 4635582bec1SHong Zhang M*/ 4645582bec1SHong Zhang 4655582bec1SHong Zhang EXTERN_C_BEGIN 4665582bec1SHong Zhang #undef __FUNCT__ 4675582bec1SHong Zhang #define __FUNCT__ "PCCreate_ML" 468dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_ML(PC pc) 4695582bec1SHong Zhang { 4705582bec1SHong Zhang PetscErrorCode ierr; 4715582bec1SHong Zhang PC_ML *pc_ml; 472776b82aeSLisandro Dalcin PetscContainer container; 4735582bec1SHong Zhang 4745582bec1SHong Zhang PetscFunctionBegin; 475573998d7SHong Zhang /* PCML is an inherited class of PCMG. Initialize pc as PCMG */ 4765582bec1SHong Zhang ierr = PCSetType(pc,PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */ 4775582bec1SHong Zhang 4785582bec1SHong Zhang /* create a supporting struct and attach it to pc */ 47938f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_ML,&pc_ml);CHKERRQ(ierr); 480776b82aeSLisandro Dalcin ierr = PetscContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr); 481776b82aeSLisandro Dalcin ierr = PetscContainerSetPointer(container,pc_ml);CHKERRQ(ierr); 482776b82aeSLisandro Dalcin ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_PC_ML);CHKERRQ(ierr); 4835582bec1SHong Zhang ierr = PetscObjectCompose((PetscObject)pc,"PC_ML",(PetscObject)container);CHKERRQ(ierr); 4845582bec1SHong Zhang 485573998d7SHong Zhang pc_ml->ml_object = 0; 486573998d7SHong Zhang pc_ml->agg_object = 0; 487573998d7SHong Zhang pc_ml->gridctx = 0; 488573998d7SHong Zhang pc_ml->PetscMLdata = 0; 489573998d7SHong Zhang pc_ml->Nlevels = -1; 490573998d7SHong Zhang pc_ml->MaxNlevels = 10; 491573998d7SHong Zhang pc_ml->MaxCoarseSize = 1; 492573998d7SHong Zhang pc_ml->CoarsenScheme = 1; /* ??? */ 493573998d7SHong Zhang pc_ml->Threshold = 0.0; 494573998d7SHong Zhang pc_ml->DampingFactor = 4.0/3.0; 495573998d7SHong Zhang pc_ml->SpectralNormScheme_Anorm = PETSC_FALSE; 496573998d7SHong Zhang pc_ml->size = 0; 497573998d7SHong Zhang 4985582bec1SHong Zhang pc_ml->PCSetUp = pc->ops->setup; 4995582bec1SHong Zhang pc_ml->PCDestroy = pc->ops->destroy; 5005582bec1SHong Zhang 5015582bec1SHong Zhang /* overwrite the pointers of PCMG by the functions of PCML */ 5025582bec1SHong Zhang pc->ops->setfromoptions = PCSetFromOptions_ML; 5035582bec1SHong Zhang pc->ops->setup = PCSetUp_ML; 5045582bec1SHong Zhang pc->ops->destroy = PCDestroy_ML; 5055582bec1SHong Zhang PetscFunctionReturn(0); 5065582bec1SHong Zhang } 5075582bec1SHong Zhang EXTERN_C_END 5085582bec1SHong Zhang 50941ca0015SHong Zhang int PetscML_getrow(ML_Operator *ML_data, int N_requested_rows, int requested_rows[], 5105582bec1SHong Zhang int allocated_space, int columns[], double values[], int row_lengths[]) 5115582bec1SHong Zhang { 5125582bec1SHong Zhang PetscErrorCode ierr; 5135582bec1SHong Zhang Mat Aloc; 5145582bec1SHong Zhang Mat_SeqAIJ *a; 5155582bec1SHong Zhang PetscInt m,i,j,k=0,row,*aj; 5165582bec1SHong Zhang PetscScalar *aa; 51741ca0015SHong Zhang FineGridCtx *ml=(FineGridCtx*)ML_Get_MyGetrowData(ML_data); 5185582bec1SHong Zhang 5195582bec1SHong Zhang Aloc = ml->Aloc; 5205582bec1SHong Zhang a = (Mat_SeqAIJ*)Aloc->data; 5215582bec1SHong Zhang ierr = MatGetSize(Aloc,&m,PETSC_NULL);CHKERRQ(ierr); 5225582bec1SHong Zhang 5235582bec1SHong Zhang for (i = 0; i<N_requested_rows; i++) { 5245582bec1SHong Zhang row = requested_rows[i]; 5255582bec1SHong Zhang row_lengths[i] = a->ilen[row]; 5265582bec1SHong Zhang if (allocated_space < k+row_lengths[i]) return(0); 5275582bec1SHong Zhang if ( (row >= 0) || (row <= (m-1)) ) { 5285582bec1SHong Zhang aj = a->j + a->i[row]; 5295582bec1SHong Zhang aa = a->a + a->i[row]; 5305582bec1SHong Zhang for (j=0; j<row_lengths[i]; j++){ 5315582bec1SHong Zhang columns[k] = aj[j]; 5325582bec1SHong Zhang values[k++] = aa[j]; 5335582bec1SHong Zhang } 5345582bec1SHong Zhang } 5355582bec1SHong Zhang } 5365582bec1SHong Zhang return(1); 5375582bec1SHong Zhang } 5385582bec1SHong Zhang 53941ca0015SHong Zhang int PetscML_matvec(ML_Operator *ML_data,int in_length,double p[],int out_length,double ap[]) 5405582bec1SHong Zhang { 5415582bec1SHong Zhang PetscErrorCode ierr; 54241ca0015SHong Zhang FineGridCtx *ml=(FineGridCtx*)ML_Get_MyMatvecData(ML_data); 5435582bec1SHong Zhang Mat A=ml->A, Aloc=ml->Aloc; 5445582bec1SHong Zhang PetscMPIInt size; 5455582bec1SHong Zhang PetscScalar *pwork=ml->pwork; 5465582bec1SHong Zhang PetscInt i; 5475582bec1SHong Zhang 5487adad957SLisandro Dalcin ierr = MPI_Comm_size(((PetscObject)A)->comm,&size);CHKERRQ(ierr); 5495582bec1SHong Zhang if (size == 1){ 55024a42b14SHong Zhang ierr = VecPlaceArray(ml->x,p);CHKERRQ(ierr); 5515582bec1SHong Zhang } else { 5525582bec1SHong Zhang for (i=0; i<in_length; i++) pwork[i] = p[i]; 5535582bec1SHong Zhang PetscML_comm(pwork,ml); 55424a42b14SHong Zhang ierr = VecPlaceArray(ml->x,pwork);CHKERRQ(ierr); 5555582bec1SHong Zhang } 55624a42b14SHong Zhang ierr = VecPlaceArray(ml->y,ap);CHKERRQ(ierr); 55724a42b14SHong Zhang ierr = MatMult(Aloc,ml->x,ml->y);CHKERRQ(ierr); 558957c941cSHong Zhang ierr = VecResetArray(ml->x);CHKERRQ(ierr); 559957c941cSHong Zhang ierr = VecResetArray(ml->y);CHKERRQ(ierr); 5605582bec1SHong Zhang return 0; 5615582bec1SHong Zhang } 5625582bec1SHong Zhang 5635582bec1SHong Zhang int PetscML_comm(double p[],void *ML_data) 5645582bec1SHong Zhang { 5655582bec1SHong Zhang PetscErrorCode ierr; 5665582bec1SHong Zhang FineGridCtx *ml=(FineGridCtx*)ML_data; 5675582bec1SHong Zhang Mat A=ml->A; 5685582bec1SHong Zhang Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data; 5695582bec1SHong Zhang PetscMPIInt size; 570a803a431SHong Zhang PetscInt i,in_length=A->rmap.n,out_length=ml->Aloc->cmap.n; 5715582bec1SHong Zhang PetscScalar *array; 5725582bec1SHong Zhang 5737adad957SLisandro Dalcin ierr = MPI_Comm_size(((PetscObject)A)->comm,&size);CHKERRQ(ierr); 5745582bec1SHong Zhang if (size == 1) return 0; 57524a42b14SHong Zhang 57624a42b14SHong Zhang ierr = VecPlaceArray(ml->y,p);CHKERRQ(ierr); 577ca9f406cSSatish Balay ierr = VecScatterBegin(a->Mvctx,ml->y,a->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 578ca9f406cSSatish Balay ierr = VecScatterEnd(a->Mvctx,ml->y,a->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 5797d15518fSHong Zhang ierr = VecResetArray(ml->y);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 } 5847d15518fSHong Zhang ierr = VecRestoreArray(a->lvec,&array);CHKERRQ(ierr); 5855582bec1SHong Zhang return 0; 5865582bec1SHong Zhang } 5875582bec1SHong Zhang #undef __FUNCT__ 5885582bec1SHong Zhang #define __FUNCT__ "MatMult_ML" 5895582bec1SHong Zhang PetscErrorCode MatMult_ML(Mat A,Vec x,Vec y) 5905582bec1SHong Zhang { 5915582bec1SHong Zhang PetscErrorCode ierr; 5925582bec1SHong Zhang Mat_MLShell *shell; 5935582bec1SHong Zhang PetscScalar *xarray,*yarray; 5945582bec1SHong Zhang PetscInt x_length,y_length; 5955582bec1SHong Zhang 5965582bec1SHong Zhang PetscFunctionBegin; 597a146b4dcSHong Zhang ierr = MatShellGetContext(A,(void **)&shell);CHKERRQ(ierr); 5985582bec1SHong Zhang ierr = VecGetArray(x,&xarray);CHKERRQ(ierr); 5995582bec1SHong Zhang ierr = VecGetArray(y,&yarray);CHKERRQ(ierr); 6005582bec1SHong Zhang x_length = shell->mlmat->invec_leng; 6015582bec1SHong Zhang y_length = shell->mlmat->outvec_leng; 6025582bec1SHong Zhang 6035582bec1SHong Zhang ML_Operator_Apply(shell->mlmat,x_length,xarray,y_length,yarray); 6045582bec1SHong Zhang 6055582bec1SHong Zhang ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr); 6065582bec1SHong Zhang ierr = VecRestoreArray(y,&yarray);CHKERRQ(ierr); 6075582bec1SHong Zhang PetscFunctionReturn(0); 6085582bec1SHong Zhang } 6095582bec1SHong Zhang /* MatMultAdd_ML - Compute y = w + A*x */ 6105582bec1SHong Zhang #undef __FUNCT__ 6115582bec1SHong Zhang #define __FUNCT__ "MatMultAdd_ML" 6125582bec1SHong Zhang PetscErrorCode MatMultAdd_ML(Mat A,Vec x,Vec w,Vec y) 6135582bec1SHong Zhang { 6145582bec1SHong Zhang PetscErrorCode ierr; 6155582bec1SHong Zhang Mat_MLShell *shell; 6165582bec1SHong Zhang PetscScalar *xarray,*yarray; 6175582bec1SHong Zhang PetscInt x_length,y_length; 6185582bec1SHong Zhang 6195582bec1SHong Zhang PetscFunctionBegin; 620a146b4dcSHong 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); 631efb30889SBarry Smith ierr = VecAXPY(y,1.0,w);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" 63975179d2cSHong Zhang PetscErrorCode MatConvert_MPIAIJ_ML(Mat A,MatType newtype,MatReuse scall,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; 646a803a431SHong Zhang PetscInt am=A->rmap.n,an=A->cmap.n,i,j,k; 6475582bec1SHong Zhang PetscInt *ci,*cj,ncols; 6485582bec1SHong Zhang 6495582bec1SHong Zhang PetscFunctionBegin; 6505582bec1SHong Zhang if (am != an) SETERRQ2(PETSC_ERR_ARG_WRONG,"A must have a square diagonal portion, am: %d != an: %d",am,an); 6515582bec1SHong Zhang 6525582bec1SHong Zhang if (scall == MAT_INITIAL_MATRIX){ 6535582bec1SHong Zhang ierr = PetscMalloc((1+am)*sizeof(PetscInt),&ci);CHKERRQ(ierr); 6545582bec1SHong Zhang ci[0] = 0; 6555582bec1SHong Zhang for (i=0; i<am; i++){ 6565582bec1SHong Zhang ci[i+1] = ci[i] + (ai[i+1] - ai[i]) + (bi[i+1] - bi[i]); 6575582bec1SHong Zhang } 6585582bec1SHong Zhang ierr = PetscMalloc((1+ci[am])*sizeof(PetscInt),&cj);CHKERRQ(ierr); 6595582bec1SHong Zhang ierr = PetscMalloc((1+ci[am])*sizeof(PetscScalar),&ca);CHKERRQ(ierr); 6605582bec1SHong Zhang 6615582bec1SHong Zhang k = 0; 6625582bec1SHong Zhang for (i=0; i<am; i++){ 6635582bec1SHong Zhang /* diagonal portion of A */ 6645582bec1SHong Zhang ncols = ai[i+1] - ai[i]; 6655582bec1SHong Zhang for (j=0; j<ncols; j++) { 6665582bec1SHong Zhang cj[k] = *aj++; 6675582bec1SHong Zhang ca[k++] = *aa++; 6685582bec1SHong Zhang } 6695582bec1SHong Zhang /* off-diagonal portion of A */ 6705582bec1SHong Zhang ncols = bi[i+1] - bi[i]; 6715582bec1SHong Zhang for (j=0; j<ncols; j++) { 6725582bec1SHong Zhang cj[k] = an + (*bj); bj++; 6735582bec1SHong Zhang ca[k++] = *ba++; 6745582bec1SHong Zhang } 6755582bec1SHong Zhang } 6765582bec1SHong Zhang if (k != ci[am]) SETERRQ2(PETSC_ERR_ARG_WRONG,"k: %d != ci[am]: %d",k,ci[am]); 6775582bec1SHong Zhang 6785582bec1SHong Zhang /* put together the new matrix */ 679a803a431SHong Zhang an = mpimat->A->cmap.n+mpimat->B->cmap.n; 6805582bec1SHong Zhang ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,an,ci,cj,ca,Aloc);CHKERRQ(ierr); 6815582bec1SHong Zhang 6825582bec1SHong Zhang /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 6835582bec1SHong Zhang /* Since these are PETSc arrays, change flags to free them as necessary. */ 6845582bec1SHong Zhang mat = (Mat_SeqAIJ*)(*Aloc)->data; 6853756052fSSatish Balay mat->free_a = PETSC_TRUE; 6863756052fSSatish Balay mat->free_ij = PETSC_TRUE; 6873756052fSSatish Balay 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; 714a146b4dcSHong 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); 718dbd8c25aSHong Zhang ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr); 7195582bec1SHong Zhang PetscFunctionReturn(0); 7205582bec1SHong Zhang } 7215582bec1SHong Zhang 7225582bec1SHong Zhang #undef __FUNCT__ 723eef31507SHong Zhang #define __FUNCT__ "MatWrapML_SeqAIJ" 724573998d7SHong Zhang PetscErrorCode MatWrapML_SeqAIJ(ML_Operator *mlmat,MatReuse reuse,Mat *newmat) 7255582bec1SHong Zhang { 726e14861a4SHong Zhang struct ML_CSR_MSRdata *matdata = (struct ML_CSR_MSRdata *)mlmat->data; 7275582bec1SHong Zhang PetscErrorCode ierr; 7283e826d7bSSatish Balay PetscInt m=mlmat->outvec_leng,n=mlmat->invec_leng,*nnz,nz_max; 729*aa85bbbfSHong Zhang PetscInt *ml_cols=matdata->columns,*ml_rowptr=matdata->rowptr,*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"); 734ab718edeSHong Zhang if (m != n){ /* ML Pmat and Rmat are in CSR format. Pass array pointers into SeqAIJ matrix */ 735573998d7SHong Zhang if (reuse){ 736573998d7SHong Zhang Mat_SeqAIJ *aij= (Mat_SeqAIJ*)(*newmat)->data; 737*aa85bbbfSHong Zhang aij->i = ml_rowptr; 738573998d7SHong Zhang aij->j = ml_cols; 739573998d7SHong Zhang aij->a = ml_vals; 740573998d7SHong Zhang } else { 741*aa85bbbfSHong Zhang /* sort ml_cols and ml_vals */ 742*aa85bbbfSHong Zhang ierr = PetscMalloc((m+1)*sizeof(PetscInt),&nnz); 743*aa85bbbfSHong Zhang for (i=0; i<m; i++) { 744*aa85bbbfSHong Zhang nnz[i] = ml_rowptr[i+1] - ml_rowptr[i]; 745*aa85bbbfSHong Zhang } 746*aa85bbbfSHong Zhang aj = ml_cols; aa = ml_vals; 747*aa85bbbfSHong Zhang for (i=0; i<m; i++){ 748*aa85bbbfSHong Zhang ierr = PetscSortIntWithScalarArray(nnz[i],aj,aa);CHKERRQ(ierr); 749*aa85bbbfSHong Zhang aj += nnz[i]; aa += nnz[i]; 750*aa85bbbfSHong Zhang } 751*aa85bbbfSHong Zhang ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,m,n,ml_rowptr,ml_cols,ml_vals,newmat);CHKERRQ(ierr); 752*aa85bbbfSHong Zhang ierr = PetscFree(nnz);CHKERRQ(ierr); 753573998d7SHong Zhang } 75424a42b14SHong Zhang PetscFunctionReturn(0); 75524a42b14SHong Zhang } 7565582bec1SHong Zhang 757e14861a4SHong Zhang /* ML Amat is in MSR format. Copy its data into SeqAIJ matrix */ 758f69a0ea3SMatthew Knepley ierr = MatCreate(PETSC_COMM_SELF,newmat);CHKERRQ(ierr); 759f69a0ea3SMatthew Knepley ierr = MatSetSizes(*newmat,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 7605582bec1SHong Zhang ierr = MatSetType(*newmat,MATSEQAIJ);CHKERRQ(ierr); 761e14861a4SHong Zhang 762573998d7SHong Zhang ierr = PetscMalloc((m+1)*sizeof(PetscInt),&nnz); 763573998d7SHong Zhang nz_max = 1; 7645582bec1SHong Zhang for (i=0; i<m; i++) { 765e14861a4SHong Zhang nnz[i] = ml_cols[i+1] - ml_cols[i] + 1; 766573998d7SHong Zhang if (nnz[i] > nz_max) nz_max += nnz[i]; 7675582bec1SHong Zhang } 7685582bec1SHong Zhang 769573998d7SHong Zhang ierr = MatSeqAIJSetPreallocation(*newmat,0,nnz);CHKERRQ(ierr); 7707c307921SBarry Smith ierr = PetscMalloc(nz_max*(sizeof(PetscInt)+sizeof(PetscScalar)),&aj);CHKERRQ(ierr); 771e14861a4SHong Zhang aa = (PetscScalar*)(aj + nz_max); 77224a42b14SHong Zhang 7735582bec1SHong Zhang for (i=0; i<m; i++){ 774e14861a4SHong Zhang k = 0; 775e14861a4SHong Zhang /* diagonal entry */ 776e14861a4SHong Zhang aj[k] = i; aa[k++] = ml_vals[i]; 777e14861a4SHong Zhang /* off diagonal entries */ 778e14861a4SHong Zhang for (j=ml_cols[i]; j<ml_cols[i+1]; j++){ 779e14861a4SHong Zhang aj[k] = ml_cols[j]; aa[k++] = ml_vals[j]; 78024a42b14SHong Zhang } 781ab718edeSHong Zhang /* sort aj and aa */ 782e14861a4SHong Zhang ierr = PetscSortIntWithScalarArray(nnz[i],aj,aa);CHKERRQ(ierr); 783e14861a4SHong Zhang ierr = MatSetValues(*newmat,1,&i,nnz[i],aj,aa,INSERT_VALUES);CHKERRQ(ierr); 7845582bec1SHong Zhang } 7855582bec1SHong Zhang ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 7865582bec1SHong Zhang ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 787573998d7SHong Zhang 788e14861a4SHong Zhang ierr = PetscFree(aj);CHKERRQ(ierr); 7893e826d7bSSatish Balay ierr = PetscFree(nnz);CHKERRQ(ierr); 7905582bec1SHong Zhang PetscFunctionReturn(0); 7915582bec1SHong Zhang } 7925582bec1SHong Zhang 7935582bec1SHong Zhang #undef __FUNCT__ 794eef31507SHong Zhang #define __FUNCT__ "MatWrapML_SHELL" 795573998d7SHong Zhang PetscErrorCode MatWrapML_SHELL(ML_Operator *mlmat,MatReuse reuse,Mat *newmat) 7965582bec1SHong Zhang { 7975582bec1SHong Zhang PetscErrorCode ierr; 7985582bec1SHong Zhang PetscInt m,n; 7995582bec1SHong Zhang ML_Comm *MLcomm; 8005582bec1SHong Zhang Mat_MLShell *shellctx; 8015582bec1SHong Zhang 8025582bec1SHong Zhang PetscFunctionBegin; 8035582bec1SHong Zhang m = mlmat->outvec_leng; 8045582bec1SHong Zhang n = mlmat->invec_leng; 8055582bec1SHong Zhang if (!m || !n){ 8065582bec1SHong Zhang newmat = PETSC_NULL; 807573998d7SHong Zhang PetscFunctionReturn(0); 808573998d7SHong Zhang } 809573998d7SHong Zhang 810573998d7SHong Zhang if (reuse){ 811573998d7SHong Zhang ierr = MatShellGetContext(*newmat,(void **)&shellctx);CHKERRQ(ierr); 812573998d7SHong Zhang shellctx->mlmat = mlmat; 813573998d7SHong Zhang PetscFunctionReturn(0); 814573998d7SHong Zhang } 815573998d7SHong Zhang 8165582bec1SHong Zhang MLcomm = mlmat->comm; 8175582bec1SHong Zhang ierr = PetscNew(Mat_MLShell,&shellctx);CHKERRQ(ierr); 8185582bec1SHong Zhang ierr = MatCreateShell(MLcomm->USR_comm,m,n,PETSC_DETERMINE,PETSC_DETERMINE,shellctx,newmat);CHKERRQ(ierr); 8193e826d7bSSatish Balay ierr = MatShellSetOperation(*newmat,MATOP_MULT,(void(*)(void))MatMult_ML);CHKERRQ(ierr); 8203e826d7bSSatish Balay ierr = MatShellSetOperation(*newmat,MATOP_MULT_ADD,(void(*)(void))MatMultAdd_ML);CHKERRQ(ierr); 8215582bec1SHong Zhang shellctx->A = *newmat; 8225582bec1SHong Zhang shellctx->mlmat = mlmat; 8235582bec1SHong Zhang ierr = VecCreate(PETSC_COMM_WORLD,&shellctx->y);CHKERRQ(ierr); 8245582bec1SHong Zhang ierr = VecSetSizes(shellctx->y,m,PETSC_DECIDE);CHKERRQ(ierr); 8255582bec1SHong Zhang ierr = VecSetFromOptions(shellctx->y);CHKERRQ(ierr); 8265582bec1SHong Zhang (*newmat)->ops->destroy = MatDestroy_ML; 8275582bec1SHong Zhang PetscFunctionReturn(0); 8285582bec1SHong Zhang } 829e14861a4SHong Zhang 830e14861a4SHong Zhang #undef __FUNCT__ 831eef31507SHong Zhang #define __FUNCT__ "MatWrapML_MPIAIJ" 832eef31507SHong Zhang PetscErrorCode MatWrapML_MPIAIJ(ML_Operator *mlmat,Mat *newmat) 833e14861a4SHong Zhang { 834ab718edeSHong Zhang struct ML_CSR_MSRdata *matdata = (struct ML_CSR_MSRdata *)mlmat->data; 835ab718edeSHong Zhang PetscInt *ml_cols=matdata->columns,*aj; 836ab718edeSHong Zhang PetscScalar *ml_vals=matdata->values,*aa; 837e14861a4SHong Zhang PetscErrorCode ierr; 838ab718edeSHong Zhang PetscInt i,j,k,*gordering; 8393e826d7bSSatish Balay PetscInt m=mlmat->outvec_leng,n,*nnzA,*nnzB,*nnz,nz_max,row; 840ab718edeSHong Zhang Mat A; 841e14861a4SHong Zhang 842e14861a4SHong Zhang PetscFunctionBegin; 843ab718edeSHong Zhang if (mlmat->getrow == NULL) SETERRQ(PETSC_ERR_ARG_NULL,"mlmat->getrow = NULL"); 844ab718edeSHong Zhang n = mlmat->invec_leng; 845e14861a4SHong Zhang if (m != n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"m %d must equal to n %d",m,n); 846e14861a4SHong Zhang 847f69a0ea3SMatthew Knepley ierr = MatCreate(mlmat->comm->USR_comm,&A);CHKERRQ(ierr); 848f69a0ea3SMatthew Knepley ierr = MatSetSizes(A,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 849ab718edeSHong Zhang ierr = MatSetType(A,MATMPIAIJ);CHKERRQ(ierr); 8503e826d7bSSatish Balay ierr = PetscMalloc3(m,PetscInt,&nnzA,m,PetscInt,&nnzB,m,PetscInt,&nnz);CHKERRQ(ierr); 8513e826d7bSSatish Balay 852e14861a4SHong Zhang nz_max = 0; 853e14861a4SHong Zhang for (i=0; i<m; i++){ 854ab718edeSHong Zhang nnz[i] = ml_cols[i+1] - ml_cols[i] + 1; 855e14861a4SHong Zhang if (nz_max < nnz[i]) nz_max = nnz[i]; 856ab718edeSHong Zhang nnzA[i] = 1; /* diag */ 857ab718edeSHong Zhang for (j=ml_cols[i]; j<ml_cols[i+1]; j++){ 858ab718edeSHong Zhang if (ml_cols[j] < m) nnzA[i]++; 859e14861a4SHong Zhang } 860e14861a4SHong Zhang nnzB[i] = nnz[i] - nnzA[i]; 861e14861a4SHong Zhang } 862ab718edeSHong Zhang ierr = MatMPIAIJSetPreallocation(A,0,nnzA,0,nnzB);CHKERRQ(ierr); 863e14861a4SHong Zhang 864ab718edeSHong Zhang /* insert mat values -- remap row and column indices */ 865ab718edeSHong Zhang nz_max++; 8667c307921SBarry Smith ierr = PetscMalloc(nz_max*(sizeof(PetscInt)+sizeof(PetscScalar)),&aj);CHKERRQ(ierr); 867ab718edeSHong Zhang aa = (PetscScalar*)(aj + nz_max); 868510c6b62SHong Zhang /* create global row numbering for a ML_Operator */ 869510c6b62SHong Zhang ML_build_global_numbering(mlmat,&gordering,"rows"); 870e14861a4SHong Zhang for (i=0; i<m; i++){ 871e14861a4SHong Zhang row = gordering[i]; 872ab718edeSHong Zhang k = 0; 873ab718edeSHong Zhang /* diagonal entry */ 874ab718edeSHong Zhang aj[k] = row; aa[k++] = ml_vals[i]; 875ab718edeSHong Zhang /* off diagonal entries */ 876ab718edeSHong Zhang for (j=ml_cols[i]; j<ml_cols[i+1]; j++){ 877ab718edeSHong Zhang aj[k] = gordering[ml_cols[j]]; aa[k++] = ml_vals[j]; 878e14861a4SHong Zhang } 879ab718edeSHong Zhang ierr = MatSetValues(A,1,&row,nnz[i],aj,aa,INSERT_VALUES);CHKERRQ(ierr); 880e14861a4SHong Zhang } 881ab718edeSHong Zhang ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 882ab718edeSHong Zhang ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 883ab718edeSHong Zhang *newmat = A; 884e14861a4SHong Zhang 8853e826d7bSSatish Balay ierr = PetscFree3(nnzA,nnzB,nnz); 886ab718edeSHong Zhang ierr = PetscFree(aj);CHKERRQ(ierr); 887e14861a4SHong Zhang PetscFunctionReturn(0); 888e14861a4SHong Zhang } 889