1dba47a55SKris Buschelman #define PETSCKSP_DLL 2ab718edeSHong Zhang 35582bec1SHong Zhang /* 42dccc152SHong Zhang Provides an interface to the ML smoothed Aggregation 57ffd031bSHong Zhang Note: Something non-obvious breaks -pc_mg_type ADDITIVE for parallel runs 67ffd031bSHong Zhang Jed Brown, see [PETSC #18321, #18449]. 75582bec1SHong Zhang */ 86356e834SBarry Smith #include "private/pcimpl.h" /*I "petscpc.h" I*/ 9*7c4f633dSBarry Smith #include "../src/ksp/pc/impls/mg/mgimpl.h" /*I "petscmg.h" I*/ 10*7c4f633dSBarry Smith #include "../src/mat/impls/aij/seq/aij.h" 11*7c4f633dSBarry Smith #include "../src/mat/impls/aij/mpi/mpiaij.h" 12cb5d8e9eSHong Zhang 135582bec1SHong Zhang #include <math.h> 142cf39c26SSatish Balay EXTERN_C_BEGIN 1568210224SSatish Balay /* HAVE_CONFIG_H flag is required by ML include files */ 1668210224SSatish Balay #if !defined(HAVE_CONFIG_H) 1768210224SSatish Balay #define HAVE_CONFIG_H 1868210224SSatish Balay #endif 195582bec1SHong Zhang #include "ml_include.h" 205582bec1SHong Zhang EXTERN_C_END 215582bec1SHong Zhang 225582bec1SHong Zhang /* The context (data structure) at each grid level */ 235582bec1SHong Zhang typedef struct { 245582bec1SHong Zhang Vec x,b,r; /* global vectors */ 255582bec1SHong Zhang Mat A,P,R; 265582bec1SHong Zhang KSP ksp; 275582bec1SHong Zhang } GridCtx; 285582bec1SHong Zhang 295582bec1SHong Zhang /* The context used to input PETSc matrix into ML at fine grid */ 305582bec1SHong Zhang typedef struct { 31573998d7SHong Zhang Mat A; /* Petsc matrix in aij format */ 32573998d7SHong Zhang Mat Aloc; /* local portion of A to be used by ML */ 3324a42b14SHong Zhang Vec x,y; 345582bec1SHong Zhang ML_Operator *mlmat; 355582bec1SHong Zhang PetscScalar *pwork; /* tmp array used by PetscML_comm() */ 365582bec1SHong Zhang } FineGridCtx; 375582bec1SHong Zhang 385582bec1SHong Zhang /* The context associates a ML matrix with a PETSc shell matrix */ 395582bec1SHong Zhang typedef struct { 405582bec1SHong Zhang Mat A; /* PETSc shell matrix associated with mlmat */ 415582bec1SHong Zhang ML_Operator *mlmat; /* ML matrix assorciated with A */ 425582bec1SHong Zhang Vec y; 435582bec1SHong Zhang } Mat_MLShell; 445582bec1SHong Zhang 455582bec1SHong Zhang /* Private context for the ML preconditioner */ 465582bec1SHong Zhang typedef struct { 475582bec1SHong Zhang ML *ml_object; 485582bec1SHong Zhang ML_Aggregate *agg_object; 495582bec1SHong Zhang GridCtx *gridctx; 505582bec1SHong Zhang FineGridCtx *PetscMLdata; 51573998d7SHong Zhang PetscInt Nlevels,MaxNlevels,MaxCoarseSize,CoarsenScheme; 525582bec1SHong Zhang PetscReal Threshold,DampingFactor; 535582bec1SHong Zhang PetscTruth SpectralNormScheme_Anorm; 54573998d7SHong Zhang PetscMPIInt size; /* size of communicator for pc->pmat */ 555582bec1SHong Zhang PetscErrorCode (*PCSetUp)(PC); 565582bec1SHong Zhang PetscErrorCode (*PCDestroy)(PC); 575582bec1SHong Zhang } PC_ML; 5841ca0015SHong Zhang 5941ca0015SHong Zhang extern int PetscML_getrow(ML_Operator *ML_data,int N_requested_rows,int requested_rows[], 605582bec1SHong Zhang int allocated_space,int columns[],double values[],int row_lengths[]); 6141ca0015SHong Zhang extern int PetscML_matvec(ML_Operator *ML_data, int in_length, double p[], int out_length,double ap[]); 625582bec1SHong Zhang extern int PetscML_comm(double x[], void *ML_data); 635582bec1SHong Zhang extern PetscErrorCode MatMult_ML(Mat,Vec,Vec); 645582bec1SHong Zhang extern PetscErrorCode MatMultAdd_ML(Mat,Vec,Vec,Vec); 6575179d2cSHong Zhang extern PetscErrorCode MatConvert_MPIAIJ_ML(Mat,MatType,MatReuse,Mat*); 665582bec1SHong Zhang extern PetscErrorCode MatDestroy_ML(Mat); 67573998d7SHong Zhang extern PetscErrorCode MatWrapML_SeqAIJ(ML_Operator*,MatReuse,Mat*); 68eef31507SHong Zhang extern PetscErrorCode MatWrapML_MPIAIJ(ML_Operator*,Mat*); 69573998d7SHong Zhang extern PetscErrorCode MatWrapML_SHELL(ML_Operator*,MatReuse,Mat*); 70573998d7SHong Zhang extern PetscErrorCode PetscContainerDestroy_PC_ML(void *); 715582bec1SHong Zhang 725582bec1SHong Zhang /* -------------------------------------------------------------------------- */ 735582bec1SHong Zhang /* 745582bec1SHong Zhang PCSetUp_ML - Prepares for the use of the ML preconditioner 755582bec1SHong Zhang by setting data structures and options. 765582bec1SHong Zhang 775582bec1SHong Zhang Input Parameter: 785582bec1SHong Zhang . pc - the preconditioner context 795582bec1SHong Zhang 805582bec1SHong Zhang Application Interface Routine: PCSetUp() 815582bec1SHong Zhang 825582bec1SHong Zhang Notes: 835582bec1SHong Zhang The interface routine PCSetUp() is not usually called directly by 845582bec1SHong Zhang the user, but instead is called by PCApply() if necessary. 855582bec1SHong Zhang */ 866ca4d86aSHong Zhang extern PetscErrorCode PCSetFromOptions_MG(PC); 875582bec1SHong Zhang #undef __FUNCT__ 885582bec1SHong Zhang #define __FUNCT__ "PCSetUp_ML" 896ca4d86aSHong Zhang PetscErrorCode PCSetUp_ML(PC pc) 905582bec1SHong Zhang { 915582bec1SHong Zhang PetscErrorCode ierr; 92eef31507SHong Zhang PetscMPIInt size; 935582bec1SHong Zhang FineGridCtx *PetscMLdata; 945582bec1SHong Zhang ML *ml_object; 955582bec1SHong Zhang ML_Aggregate *agg_object; 965582bec1SHong Zhang ML_Operator *mlmat; 97ac346b81SHong Zhang PetscInt nlocal_allcols,Nlevels,mllevel,level,level1,m,fine_level; 985582bec1SHong Zhang Mat A,Aloc; 995582bec1SHong Zhang GridCtx *gridctx; 1005582bec1SHong Zhang PC_ML *pc_ml=PETSC_NULL; 101776b82aeSLisandro Dalcin PetscContainer container; 102573998d7SHong Zhang MatReuse reuse = MAT_INITIAL_MATRIX; 1035582bec1SHong Zhang 1045582bec1SHong Zhang PetscFunctionBegin; 1055582bec1SHong Zhang ierr = PetscObjectQuery((PetscObject)pc,"PC_ML",(PetscObject *)&container);CHKERRQ(ierr); 1065582bec1SHong Zhang if (container) { 107776b82aeSLisandro Dalcin ierr = PetscContainerGetPointer(container,(void **)&pc_ml);CHKERRQ(ierr); 1085582bec1SHong Zhang } else { 1095582bec1SHong Zhang SETERRQ(PETSC_ERR_ARG_NULL,"Container does not exit"); 1105582bec1SHong Zhang } 1115582bec1SHong Zhang 112573998d7SHong Zhang if (pc->setupcalled){ 113573998d7SHong Zhang if (pc->flag == SAME_NONZERO_PATTERN){ 114573998d7SHong Zhang reuse = MAT_REUSE_MATRIX; 115573998d7SHong Zhang PetscMLdata = pc_ml->PetscMLdata; 116573998d7SHong Zhang gridctx = pc_ml->gridctx; 117573998d7SHong Zhang /* ML objects cannot be reused */ 118573998d7SHong Zhang ML_Destroy(&pc_ml->ml_object); 119573998d7SHong Zhang ML_Aggregate_Destroy(&pc_ml->agg_object); 120573998d7SHong Zhang } else { 121573998d7SHong Zhang PC_ML *pc_ml_new = PETSC_NULL; 122573998d7SHong Zhang PetscContainer container_new; 12338f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_ML,&pc_ml_new);CHKERRQ(ierr); 124573998d7SHong Zhang ierr = PetscContainerCreate(PETSC_COMM_SELF,&container_new);CHKERRQ(ierr); 125573998d7SHong Zhang ierr = PetscContainerSetPointer(container_new,pc_ml_new);CHKERRQ(ierr); 126573998d7SHong Zhang ierr = PetscContainerSetUserDestroy(container_new,PetscContainerDestroy_PC_ML);CHKERRQ(ierr); 127573998d7SHong Zhang ierr = PetscObjectCompose((PetscObject)pc,"PC_ML",(PetscObject)container_new);CHKERRQ(ierr); 128573998d7SHong Zhang 129573998d7SHong Zhang ierr = PetscMemcpy(pc_ml_new,pc_ml,sizeof(PC_ML));CHKERRQ(ierr); 130573998d7SHong Zhang ierr = PetscContainerDestroy(container);CHKERRQ(ierr); 131573998d7SHong Zhang pc_ml = pc_ml_new; 132573998d7SHong Zhang } 133573998d7SHong Zhang } 134573998d7SHong Zhang 1355582bec1SHong Zhang /* setup special features of PCML */ 1365582bec1SHong Zhang /*--------------------------------*/ 1375582bec1SHong Zhang /* covert A to Aloc to be used by ML at fine grid */ 1385582bec1SHong Zhang A = pc->pmat; 1397adad957SLisandro Dalcin ierr = MPI_Comm_size(((PetscObject)A)->comm,&size);CHKERRQ(ierr); 1405582bec1SHong Zhang pc_ml->size = size; 1415582bec1SHong Zhang if (size > 1){ 142573998d7SHong Zhang if (reuse) Aloc = PetscMLdata->Aloc; 143573998d7SHong Zhang ierr = MatConvert_MPIAIJ_ML(A,PETSC_NULL,reuse,&Aloc);CHKERRQ(ierr); 1445582bec1SHong Zhang } else { 1455582bec1SHong Zhang Aloc = A; 1465582bec1SHong Zhang } 1475582bec1SHong Zhang 1485582bec1SHong Zhang /* create and initialize struct 'PetscMLdata' */ 149573998d7SHong Zhang if (!reuse){ 15038f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,FineGridCtx,&PetscMLdata);CHKERRQ(ierr); 1515582bec1SHong Zhang pc_ml->PetscMLdata = PetscMLdata; 152d0f46423SBarry Smith ierr = PetscMalloc((Aloc->cmap->n+1)*sizeof(PetscScalar),&PetscMLdata->pwork);CHKERRQ(ierr); 1535582bec1SHong Zhang 15424a42b14SHong Zhang ierr = VecCreate(PETSC_COMM_SELF,&PetscMLdata->x);CHKERRQ(ierr); 155d0f46423SBarry Smith ierr = VecSetSizes(PetscMLdata->x,Aloc->cmap->n,Aloc->cmap->n);CHKERRQ(ierr); 15624a42b14SHong Zhang ierr = VecSetType(PetscMLdata->x,VECSEQ);CHKERRQ(ierr); 15724a42b14SHong Zhang 15824a42b14SHong Zhang ierr = VecCreate(PETSC_COMM_SELF,&PetscMLdata->y);CHKERRQ(ierr); 159d0f46423SBarry Smith ierr = VecSetSizes(PetscMLdata->y,A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr); 16024a42b14SHong Zhang ierr = VecSetType(PetscMLdata->y,VECSEQ);CHKERRQ(ierr); 161573998d7SHong Zhang } 162573998d7SHong Zhang PetscMLdata->A = A; 163573998d7SHong Zhang PetscMLdata->Aloc = Aloc; 16424a42b14SHong Zhang 1655582bec1SHong Zhang /* create ML discretization matrix at fine grid */ 16645cf47abSHong Zhang /* ML requires input of fine-grid matrix. It determines nlevels. */ 1675582bec1SHong Zhang ierr = MatGetSize(Aloc,&m,&nlocal_allcols);CHKERRQ(ierr); 1685582bec1SHong Zhang ML_Create(&ml_object,pc_ml->MaxNlevels); 169573998d7SHong Zhang pc_ml->ml_object = ml_object; 1705582bec1SHong Zhang ML_Init_Amatrix(ml_object,0,m,m,PetscMLdata); 1715582bec1SHong Zhang ML_Set_Amatrix_Getrow(ml_object,0,PetscML_getrow,PetscML_comm,nlocal_allcols); 1725582bec1SHong Zhang ML_Set_Amatrix_Matvec(ml_object,0,PetscML_matvec); 1735582bec1SHong Zhang 1745582bec1SHong Zhang /* aggregation */ 1755582bec1SHong Zhang ML_Aggregate_Create(&agg_object); 176573998d7SHong Zhang pc_ml->agg_object = agg_object; 177573998d7SHong Zhang 1785582bec1SHong Zhang ML_Aggregate_Set_MaxCoarseSize(agg_object,pc_ml->MaxCoarseSize); 1795582bec1SHong Zhang /* set options */ 1805582bec1SHong Zhang switch (pc_ml->CoarsenScheme) { 1815582bec1SHong Zhang case 1: 1825582bec1SHong Zhang ML_Aggregate_Set_CoarsenScheme_Coupled(agg_object);break; 1835582bec1SHong Zhang case 2: 1845582bec1SHong Zhang ML_Aggregate_Set_CoarsenScheme_MIS(agg_object);break; 1855582bec1SHong Zhang case 3: 1865582bec1SHong Zhang ML_Aggregate_Set_CoarsenScheme_METIS(agg_object);break; 1875582bec1SHong Zhang } 1885582bec1SHong Zhang ML_Aggregate_Set_Threshold(agg_object,pc_ml->Threshold); 1895582bec1SHong Zhang ML_Aggregate_Set_DampingFactor(agg_object,pc_ml->DampingFactor); 1905582bec1SHong Zhang if (pc_ml->SpectralNormScheme_Anorm){ 1917ffd031bSHong Zhang ML_Set_SpectralNormScheme_Anorm(ml_object); 1925582bec1SHong Zhang } 1935582bec1SHong Zhang 1945582bec1SHong Zhang Nlevels = ML_Gen_MGHierarchy_UsingAggregation(ml_object,0,ML_INCREASING,agg_object); 1955582bec1SHong Zhang if (Nlevels<=0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Nlevels %d must > 0",Nlevels); 196573998d7SHong 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); 197573998d7SHong Zhang pc_ml->Nlevels = Nlevels; 198aa85bbbfSHong Zhang fine_level = Nlevels - 1; 199573998d7SHong Zhang if (!pc->setupcalled){ 20097177400SBarry Smith ierr = PCMGSetLevels(pc,Nlevels,PETSC_NULL);CHKERRQ(ierr); 201aa85bbbfSHong Zhang /* set default smoothers */ 202aa85bbbfSHong Zhang KSP smoother; 203aa85bbbfSHong Zhang PC subpc; 204aa85bbbfSHong Zhang for (level=1; level<=fine_level; level++){ 205aa85bbbfSHong Zhang if (size == 1){ 206aa85bbbfSHong Zhang ierr = PCMGGetSmoother(pc,level,&smoother);CHKERRQ(ierr); 207aa85bbbfSHong Zhang ierr = KSPSetType(smoother,KSPRICHARDSON);CHKERRQ(ierr); 208aa85bbbfSHong Zhang ierr = KSPGetPC(smoother,&subpc);CHKERRQ(ierr); 209aa85bbbfSHong Zhang ierr = PCSetType(subpc,PCSOR);CHKERRQ(ierr); 210aa85bbbfSHong Zhang } else { 211aa85bbbfSHong Zhang ierr = PCMGGetSmoother(pc,level,&smoother);CHKERRQ(ierr); 212aa85bbbfSHong Zhang ierr = KSPSetType(smoother,KSPRICHARDSON);CHKERRQ(ierr); 213aa85bbbfSHong Zhang ierr = KSPGetPC(smoother,&subpc);CHKERRQ(ierr); 214aa85bbbfSHong Zhang ierr = PCSetType(subpc,PCSOR);CHKERRQ(ierr); 215aa85bbbfSHong Zhang } 216aa85bbbfSHong Zhang } 21797177400SBarry Smith ierr = PCSetFromOptions_MG(pc);CHKERRQ(ierr); /* should be called in PCSetFromOptions_ML(), but cannot be called prior to PCMGSetLevels() */ 218573998d7SHong Zhang } 2195582bec1SHong Zhang 220573998d7SHong Zhang if (!reuse){ 2215582bec1SHong Zhang ierr = PetscMalloc(Nlevels*sizeof(GridCtx),&gridctx);CHKERRQ(ierr); 2225582bec1SHong Zhang pc_ml->gridctx = gridctx; 223573998d7SHong Zhang } 2245582bec1SHong Zhang 2255582bec1SHong Zhang /* wrap ML matrices by PETSc shell matrices at coarsened grids. 2265582bec1SHong Zhang Level 0 is the finest grid for ML, but coarsest for PETSc! */ 227e14861a4SHong Zhang gridctx[fine_level].A = A; 228573998d7SHong Zhang 229e14861a4SHong Zhang level = fine_level - 1; 230ab718edeSHong Zhang if (size == 1){ /* convert ML P, R and A into seqaij format */ 2315582bec1SHong Zhang for (mllevel=1; mllevel<Nlevels; mllevel++){ 232e14861a4SHong Zhang mlmat = &(ml_object->Pmat[mllevel]); 233573998d7SHong Zhang ierr = MatWrapML_SeqAIJ(mlmat,reuse,&gridctx[level].P);CHKERRQ(ierr); 234e14861a4SHong Zhang mlmat = &(ml_object->Rmat[mllevel-1]); 235573998d7SHong Zhang ierr = MatWrapML_SeqAIJ(mlmat,reuse,&gridctx[level].R);CHKERRQ(ierr); 236573998d7SHong Zhang 237573998d7SHong Zhang mlmat = &(ml_object->Amat[mllevel]); 238573998d7SHong Zhang if (reuse){ 239573998d7SHong Zhang /* ML matrix A changes sparse pattern although PETSc A doesn't, thus gridctx[level].A must be recreated! */ 240573998d7SHong Zhang ierr = MatDestroy(gridctx[level].A);CHKERRQ(ierr); 241573998d7SHong Zhang } 242573998d7SHong Zhang ierr = MatWrapML_SeqAIJ(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].A);CHKERRQ(ierr); 2435582bec1SHong Zhang level--; 2445582bec1SHong Zhang } 245ab718edeSHong Zhang } else { /* convert ML P and R into shell format, ML A into mpiaij format */ 2465582bec1SHong Zhang for (mllevel=1; mllevel<Nlevels; mllevel++){ 2475582bec1SHong Zhang mlmat = &(ml_object->Pmat[mllevel]); 248573998d7SHong Zhang ierr = MatWrapML_SHELL(mlmat,reuse,&gridctx[level].P);CHKERRQ(ierr); 249ab718edeSHong Zhang mlmat = &(ml_object->Rmat[mllevel-1]); 250573998d7SHong Zhang ierr = MatWrapML_SHELL(mlmat,reuse,&gridctx[level].R);CHKERRQ(ierr); 251573998d7SHong Zhang 2525582bec1SHong Zhang mlmat = &(ml_object->Amat[mllevel]); 253573998d7SHong Zhang if (reuse){ 254573998d7SHong Zhang ierr = MatDestroy(gridctx[level].A);CHKERRQ(ierr); 255573998d7SHong Zhang } 256eef31507SHong Zhang ierr = MatWrapML_MPIAIJ(mlmat,&gridctx[level].A);CHKERRQ(ierr); 2575582bec1SHong Zhang level--; 2585582bec1SHong Zhang } 2595582bec1SHong Zhang } 2605582bec1SHong Zhang 261573998d7SHong Zhang /* create vectors and ksp at all levels */ 262573998d7SHong Zhang if (!reuse){ 263ac346b81SHong Zhang for (level=0; level<fine_level; level++){ 264573998d7SHong Zhang level1 = level + 1; 265e64afeacSLisandro Dalcin ierr = VecCreate(((PetscObject)gridctx[level].A)->comm,&gridctx[level].x);CHKERRQ(ierr); 266d0f46423SBarry Smith ierr = VecSetSizes(gridctx[level].x,gridctx[level].A->cmap->n,PETSC_DECIDE);CHKERRQ(ierr); 2675582bec1SHong Zhang ierr = VecSetType(gridctx[level].x,VECMPI);CHKERRQ(ierr); 26897177400SBarry Smith ierr = PCMGSetX(pc,level,gridctx[level].x);CHKERRQ(ierr); 2695582bec1SHong Zhang 270e64afeacSLisandro Dalcin ierr = VecCreate(((PetscObject)gridctx[level].A)->comm,&gridctx[level].b);CHKERRQ(ierr); 271d0f46423SBarry Smith ierr = VecSetSizes(gridctx[level].b,gridctx[level].A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr); 2725582bec1SHong Zhang ierr = VecSetType(gridctx[level].b,VECMPI);CHKERRQ(ierr); 27397177400SBarry Smith ierr = PCMGSetRhs(pc,level,gridctx[level].b);CHKERRQ(ierr); 274ac346b81SHong Zhang 275e64afeacSLisandro Dalcin ierr = VecCreate(((PetscObject)gridctx[level1].A)->comm,&gridctx[level1].r);CHKERRQ(ierr); 276d0f46423SBarry Smith ierr = VecSetSizes(gridctx[level1].r,gridctx[level1].A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr); 277ac346b81SHong Zhang ierr = VecSetType(gridctx[level1].r,VECMPI);CHKERRQ(ierr); 27897177400SBarry Smith ierr = PCMGSetR(pc,level1,gridctx[level1].r);CHKERRQ(ierr); 279ac346b81SHong Zhang 2805582bec1SHong Zhang if (level == 0){ 28197177400SBarry Smith ierr = PCMGGetCoarseSolve(pc,&gridctx[level].ksp);CHKERRQ(ierr); 2825582bec1SHong Zhang } else { 28397177400SBarry Smith ierr = PCMGGetSmoother(pc,level,&gridctx[level].ksp);CHKERRQ(ierr); 284573998d7SHong Zhang } 285573998d7SHong Zhang } 286573998d7SHong Zhang ierr = PCMGGetSmoother(pc,fine_level,&gridctx[fine_level].ksp);CHKERRQ(ierr); 287573998d7SHong Zhang } 288573998d7SHong Zhang 289573998d7SHong Zhang /* create coarse level and the interpolation between the levels */ 290573998d7SHong Zhang for (level=0; level<fine_level; level++){ 291573998d7SHong Zhang level1 = level + 1; 292aea2a34eSBarry Smith ierr = PCMGSetInterpolation(pc,level1,gridctx[level].P);CHKERRQ(ierr); 293573998d7SHong Zhang ierr = PCMGSetRestriction(pc,level1,gridctx[level].R);CHKERRQ(ierr); 294573998d7SHong Zhang if (level > 0){ 29597177400SBarry Smith ierr = PCMGSetResidual(pc,level,PCMGDefaultResidual,gridctx[level].A);CHKERRQ(ierr); 2965582bec1SHong Zhang } 2975582bec1SHong Zhang ierr = KSPSetOperators(gridctx[level].ksp,gridctx[level].A,gridctx[level].A,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 2985582bec1SHong Zhang } 29997177400SBarry Smith ierr = PCMGSetResidual(pc,fine_level,PCMGDefaultResidual,gridctx[fine_level].A);CHKERRQ(ierr); 300ac346b81SHong Zhang ierr = KSPSetOperators(gridctx[fine_level].ksp,gridctx[level].A,gridctx[fine_level].A,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 3015582bec1SHong Zhang 3025582bec1SHong Zhang /* now call PCSetUp_MG() */ 303573998d7SHong Zhang /*-------------------------------*/ 3045582bec1SHong Zhang ierr = (*pc_ml->PCSetUp)(pc);CHKERRQ(ierr); 3055582bec1SHong Zhang PetscFunctionReturn(0); 3065582bec1SHong Zhang } 3075582bec1SHong Zhang 3085582bec1SHong Zhang #undef __FUNCT__ 309776b82aeSLisandro Dalcin #define __FUNCT__ "PetscContainerDestroy_PC_ML" 310776b82aeSLisandro Dalcin PetscErrorCode PetscContainerDestroy_PC_ML(void *ptr) 3115582bec1SHong Zhang { 3125582bec1SHong Zhang PetscErrorCode ierr; 3135582bec1SHong Zhang PC_ML *pc_ml = (PC_ML*)ptr; 314573998d7SHong Zhang PetscInt level,fine_level=pc_ml->Nlevels-1; 3155582bec1SHong Zhang 3165582bec1SHong Zhang PetscFunctionBegin; 3175582bec1SHong Zhang ML_Aggregate_Destroy(&pc_ml->agg_object); 3185582bec1SHong Zhang ML_Destroy(&pc_ml->ml_object); 3195582bec1SHong Zhang 32038f2d2fdSLisandro Dalcin if (pc_ml->PetscMLdata) { 3215582bec1SHong Zhang ierr = PetscFree(pc_ml->PetscMLdata->pwork);CHKERRQ(ierr); 32238f2d2fdSLisandro Dalcin if (pc_ml->size > 1) {ierr = MatDestroy(pc_ml->PetscMLdata->Aloc);CHKERRQ(ierr);} 323ac346b81SHong Zhang if (pc_ml->PetscMLdata->x){ierr = VecDestroy(pc_ml->PetscMLdata->x);CHKERRQ(ierr);} 324ac346b81SHong Zhang if (pc_ml->PetscMLdata->y){ierr = VecDestroy(pc_ml->PetscMLdata->y);CHKERRQ(ierr);} 32538f2d2fdSLisandro Dalcin } 3265582bec1SHong Zhang ierr = PetscFree(pc_ml->PetscMLdata);CHKERRQ(ierr); 3275582bec1SHong Zhang 328573998d7SHong Zhang for (level=0; level<fine_level; level++){ 329ac346b81SHong Zhang if (pc_ml->gridctx[level].A){ierr = MatDestroy(pc_ml->gridctx[level].A);CHKERRQ(ierr);} 330ac346b81SHong Zhang if (pc_ml->gridctx[level].P){ierr = MatDestroy(pc_ml->gridctx[level].P);CHKERRQ(ierr);} 331ac346b81SHong Zhang if (pc_ml->gridctx[level].R){ierr = MatDestroy(pc_ml->gridctx[level].R);CHKERRQ(ierr);} 332ac346b81SHong Zhang if (pc_ml->gridctx[level].x){ierr = VecDestroy(pc_ml->gridctx[level].x);CHKERRQ(ierr);} 333ac346b81SHong Zhang if (pc_ml->gridctx[level].b){ierr = VecDestroy(pc_ml->gridctx[level].b);CHKERRQ(ierr);} 334ac346b81SHong Zhang if (pc_ml->gridctx[level+1].r){ierr = VecDestroy(pc_ml->gridctx[level+1].r);CHKERRQ(ierr);} 3355582bec1SHong Zhang } 3365582bec1SHong Zhang ierr = PetscFree(pc_ml->gridctx);CHKERRQ(ierr); 3375582bec1SHong Zhang ierr = PetscFree(pc_ml);CHKERRQ(ierr); 3385582bec1SHong Zhang PetscFunctionReturn(0); 3395582bec1SHong Zhang } 3405582bec1SHong Zhang /* -------------------------------------------------------------------------- */ 3415582bec1SHong Zhang /* 3425582bec1SHong Zhang PCDestroy_ML - Destroys the private context for the ML preconditioner 3435582bec1SHong Zhang that was created with PCCreate_ML(). 3445582bec1SHong Zhang 3455582bec1SHong Zhang Input Parameter: 3465582bec1SHong Zhang . pc - the preconditioner context 3475582bec1SHong Zhang 3485582bec1SHong Zhang Application Interface Routine: PCDestroy() 3495582bec1SHong Zhang */ 3505582bec1SHong Zhang #undef __FUNCT__ 3515582bec1SHong Zhang #define __FUNCT__ "PCDestroy_ML" 3526ca4d86aSHong Zhang PetscErrorCode PCDestroy_ML(PC pc) 3535582bec1SHong Zhang { 3545582bec1SHong Zhang PetscErrorCode ierr; 3555582bec1SHong Zhang PC_ML *pc_ml=PETSC_NULL; 356776b82aeSLisandro Dalcin PetscContainer container; 3575582bec1SHong Zhang 3585582bec1SHong Zhang PetscFunctionBegin; 3595582bec1SHong Zhang ierr = PetscObjectQuery((PetscObject)pc,"PC_ML",(PetscObject *)&container);CHKERRQ(ierr); 3605582bec1SHong Zhang if (container) { 361776b82aeSLisandro Dalcin ierr = PetscContainerGetPointer(container,(void **)&pc_ml);CHKERRQ(ierr); 3625582bec1SHong Zhang pc->ops->destroy = pc_ml->PCDestroy; 3635582bec1SHong Zhang } else { 3645582bec1SHong Zhang SETERRQ(PETSC_ERR_ARG_NULL,"Container does not exit"); 3655582bec1SHong Zhang } 3669cb0ec93SHong Zhang /* detach pc and PC_ML and dereference container */ 36738f2d2fdSLisandro Dalcin ierr = PetscContainerDestroy(container);CHKERRQ(ierr); 3685582bec1SHong Zhang ierr = PetscObjectCompose((PetscObject)pc,"PC_ML",0);CHKERRQ(ierr); 36938f2d2fdSLisandro Dalcin if (pc->ops->destroy) { 3705582bec1SHong Zhang ierr = (*pc->ops->destroy)(pc);CHKERRQ(ierr); 37138f2d2fdSLisandro Dalcin } 3725582bec1SHong Zhang PetscFunctionReturn(0); 3735582bec1SHong Zhang } 3745582bec1SHong Zhang 3755582bec1SHong Zhang #undef __FUNCT__ 3765582bec1SHong Zhang #define __FUNCT__ "PCSetFromOptions_ML" 3776ca4d86aSHong Zhang PetscErrorCode PCSetFromOptions_ML(PC pc) 3785582bec1SHong Zhang { 3795582bec1SHong Zhang PetscErrorCode ierr; 38038f2d2fdSLisandro Dalcin PetscInt indx,m,PrintLevel; 3815582bec1SHong Zhang PetscTruth flg; 3825582bec1SHong Zhang const char *scheme[] = {"Uncoupled","Coupled","MIS","METIS"}; 3835582bec1SHong Zhang PC_ML *pc_ml=PETSC_NULL; 384776b82aeSLisandro Dalcin PetscContainer container; 3859dcbbd2bSBarry Smith PCMGType mgtype; 3865582bec1SHong Zhang 3875582bec1SHong Zhang PetscFunctionBegin; 3885582bec1SHong Zhang ierr = PetscObjectQuery((PetscObject)pc,"PC_ML",(PetscObject *)&container);CHKERRQ(ierr); 3895582bec1SHong Zhang if (container) { 390776b82aeSLisandro Dalcin ierr = PetscContainerGetPointer(container,(void **)&pc_ml);CHKERRQ(ierr); 3915582bec1SHong Zhang } else { 3925582bec1SHong Zhang SETERRQ(PETSC_ERR_ARG_NULL,"Container does not exit"); 3935582bec1SHong Zhang } 3946ca4d86aSHong Zhang 3955582bec1SHong Zhang /* inherited MG options */ 3966ca4d86aSHong Zhang ierr = PetscOptionsHead("Multigrid options(inherited)");CHKERRQ(ierr); 3975582bec1SHong Zhang ierr = PetscOptionsInt("-pc_mg_cycles","1 for V cycle, 2 for W-cycle","MGSetCycles",1,&m,&flg);CHKERRQ(ierr); 3985582bec1SHong Zhang ierr = PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","MGSetNumberSmoothUp",1,&m,&flg);CHKERRQ(ierr); 3995582bec1SHong Zhang ierr = PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","MGSetNumberSmoothDown",1,&m,&flg);CHKERRQ(ierr); 4009dcbbd2bSBarry Smith ierr = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)PC_MG_MULTIPLICATIVE,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr); 4015582bec1SHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 4025582bec1SHong Zhang 4035582bec1SHong Zhang /* ML options */ 4045582bec1SHong Zhang ierr = PetscOptionsHead("ML options");CHKERRQ(ierr); 4055582bec1SHong Zhang /* set defaults */ 4065582bec1SHong Zhang PrintLevel = 0; 4075582bec1SHong Zhang indx = 0; 4085582bec1SHong Zhang ierr = PetscOptionsInt("-pc_ml_PrintLevel","Print level","ML_Set_PrintLevel",PrintLevel,&PrintLevel,PETSC_NULL);CHKERRQ(ierr); 4095582bec1SHong Zhang ML_Set_PrintLevel(PrintLevel); 410573998d7SHong Zhang ierr = PetscOptionsInt("-pc_ml_maxNlevels","Maximum number of levels","None",pc_ml->MaxNlevels,&pc_ml->MaxNlevels,PETSC_NULL);CHKERRQ(ierr); 411573998d7SHong Zhang ierr = PetscOptionsInt("-pc_ml_maxCoarseSize","Maximum coarsest mesh size","ML_Aggregate_Set_MaxCoarseSize",pc_ml->MaxCoarseSize,&pc_ml->MaxCoarseSize,PETSC_NULL);CHKERRQ(ierr); 412573998d7SHong Zhang ierr = PetscOptionsEList("-pc_ml_CoarsenScheme","Aggregate Coarsen Scheme","ML_Aggregate_Set_CoarsenScheme_*",scheme,4,scheme[0],&indx,PETSC_NULL);CHKERRQ(ierr); /* ??? */ 4135582bec1SHong Zhang pc_ml->CoarsenScheme = indx; 4145582bec1SHong Zhang 415573998d7SHong Zhang ierr = PetscOptionsReal("-pc_ml_DampingFactor","P damping factor","ML_Aggregate_Set_DampingFactor",pc_ml->DampingFactor,&pc_ml->DampingFactor,PETSC_NULL);CHKERRQ(ierr); 4165582bec1SHong Zhang 417573998d7SHong Zhang ierr = PetscOptionsReal("-pc_ml_Threshold","Smoother drop tol","ML_Aggregate_Set_Threshold",pc_ml->Threshold,&pc_ml->Threshold,PETSC_NULL);CHKERRQ(ierr); 4185582bec1SHong Zhang 41940597110SHong Zhang ierr = PetscOptionsTruth("-pc_ml_SpectralNormScheme_Anorm","Method used for estimating spectral radius","ML_Set_SpectralNormScheme_Anorm",pc_ml->SpectralNormScheme_Anorm,&pc_ml->SpectralNormScheme_Anorm,PETSC_NULL); 4205582bec1SHong Zhang 4215582bec1SHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 4225582bec1SHong Zhang PetscFunctionReturn(0); 4235582bec1SHong Zhang } 4245582bec1SHong Zhang 4255582bec1SHong Zhang /* -------------------------------------------------------------------------- */ 4265582bec1SHong Zhang /* 4275582bec1SHong Zhang PCCreate_ML - Creates a ML preconditioner context, PC_ML, 4285582bec1SHong Zhang and sets this as the private data within the generic preconditioning 4295582bec1SHong Zhang context, PC, that was created within PCCreate(). 4305582bec1SHong Zhang 4315582bec1SHong Zhang Input Parameter: 4325582bec1SHong Zhang . pc - the preconditioner context 4335582bec1SHong Zhang 4345582bec1SHong Zhang Application Interface Routine: PCCreate() 4355582bec1SHong Zhang */ 4365582bec1SHong Zhang 4375582bec1SHong Zhang /*MC 4381e5ab15bSHong Zhang PCML - Use algebraic multigrid preconditioning. This preconditioner requires you provide 4395582bec1SHong Zhang fine grid discretization matrix. The coarser grid matrices and restriction/interpolation 4406ca4d86aSHong Zhang operators are computed by ML, with the matrices coverted to PETSc matrices in aij format 4416ca4d86aSHong Zhang and the restriction/interpolation operators wrapped as PETSc shell matrices. 4425582bec1SHong Zhang 4436ca4d86aSHong Zhang Options Database Key: 4446ca4d86aSHong Zhang Multigrid options(inherited) 4456ca4d86aSHong Zhang + -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (MGSetCycles) 4466ca4d86aSHong Zhang . -pc_mg_smoothup <1>: Number of post-smoothing steps (MGSetNumberSmoothUp) 4476ca4d86aSHong Zhang . -pc_mg_smoothdown <1>: Number of pre-smoothing steps (MGSetNumberSmoothDown) 448f41ab451SVictor Eijkhout - -pc_mg_type <multiplicative>: (one of) additive multiplicative full cascade kascade 4496ca4d86aSHong Zhang 45051f519a2SBarry Smith ML options: 4516ca4d86aSHong Zhang + -pc_ml_PrintLevel <0>: Print level (ML_Set_PrintLevel) 4526ca4d86aSHong Zhang . -pc_ml_maxNlevels <10>: Maximum number of levels (None) 4536ca4d86aSHong Zhang . -pc_ml_maxCoarseSize <1>: Maximum coarsest mesh size (ML_Aggregate_Set_MaxCoarseSize) 454f41ab451SVictor Eijkhout . -pc_ml_CoarsenScheme <Uncoupled>: (one of) Uncoupled Coupled MIS METIS 4556ca4d86aSHong Zhang . -pc_ml_DampingFactor <1.33333>: P damping factor (ML_Aggregate_Set_DampingFactor) 4566ca4d86aSHong Zhang . -pc_ml_Threshold <0>: Smoother drop tol (ML_Aggregate_Set_Threshold) 4577ffd031bSHong Zhang - -pc_ml_SpectralNormScheme_Anorm <false>: Method used for estimating spectral radius (ML_Set_SpectralNormScheme_Anorm) 4585582bec1SHong Zhang 4595582bec1SHong Zhang Level: intermediate 4605582bec1SHong Zhang 4615582bec1SHong Zhang Concepts: multigrid 4625582bec1SHong Zhang 4635582bec1SHong Zhang .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, 46497177400SBarry Smith PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), MPSetCycles(), PCMGSetNumberSmoothDown(), 46597177400SBarry Smith PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(), 46697177400SBarry Smith PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(), 46797177400SBarry Smith PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR() 4685582bec1SHong Zhang M*/ 4695582bec1SHong Zhang 4705582bec1SHong Zhang EXTERN_C_BEGIN 4715582bec1SHong Zhang #undef __FUNCT__ 4725582bec1SHong Zhang #define __FUNCT__ "PCCreate_ML" 473dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_ML(PC pc) 4745582bec1SHong Zhang { 4755582bec1SHong Zhang PetscErrorCode ierr; 4765582bec1SHong Zhang PC_ML *pc_ml; 477776b82aeSLisandro Dalcin PetscContainer container; 4785582bec1SHong Zhang 4795582bec1SHong Zhang PetscFunctionBegin; 480573998d7SHong Zhang /* PCML is an inherited class of PCMG. Initialize pc as PCMG */ 4815582bec1SHong Zhang ierr = PCSetType(pc,PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */ 4825582bec1SHong Zhang 4835582bec1SHong Zhang /* create a supporting struct and attach it to pc */ 48438f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_ML,&pc_ml);CHKERRQ(ierr); 485776b82aeSLisandro Dalcin ierr = PetscContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr); 486776b82aeSLisandro Dalcin ierr = PetscContainerSetPointer(container,pc_ml);CHKERRQ(ierr); 487776b82aeSLisandro Dalcin ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_PC_ML);CHKERRQ(ierr); 4885582bec1SHong Zhang ierr = PetscObjectCompose((PetscObject)pc,"PC_ML",(PetscObject)container);CHKERRQ(ierr); 4895582bec1SHong Zhang 490573998d7SHong Zhang pc_ml->ml_object = 0; 491573998d7SHong Zhang pc_ml->agg_object = 0; 492573998d7SHong Zhang pc_ml->gridctx = 0; 493573998d7SHong Zhang pc_ml->PetscMLdata = 0; 494573998d7SHong Zhang pc_ml->Nlevels = -1; 495573998d7SHong Zhang pc_ml->MaxNlevels = 10; 496573998d7SHong Zhang pc_ml->MaxCoarseSize = 1; 497573998d7SHong Zhang pc_ml->CoarsenScheme = 1; /* ??? */ 498573998d7SHong Zhang pc_ml->Threshold = 0.0; 499573998d7SHong Zhang pc_ml->DampingFactor = 4.0/3.0; 500573998d7SHong Zhang pc_ml->SpectralNormScheme_Anorm = PETSC_FALSE; 501573998d7SHong Zhang pc_ml->size = 0; 502573998d7SHong Zhang 5035582bec1SHong Zhang pc_ml->PCSetUp = pc->ops->setup; 5045582bec1SHong Zhang pc_ml->PCDestroy = pc->ops->destroy; 5055582bec1SHong Zhang 5065582bec1SHong Zhang /* overwrite the pointers of PCMG by the functions of PCML */ 5075582bec1SHong Zhang pc->ops->setfromoptions = PCSetFromOptions_ML; 5085582bec1SHong Zhang pc->ops->setup = PCSetUp_ML; 5095582bec1SHong Zhang pc->ops->destroy = PCDestroy_ML; 5105582bec1SHong Zhang PetscFunctionReturn(0); 5115582bec1SHong Zhang } 5125582bec1SHong Zhang EXTERN_C_END 5135582bec1SHong Zhang 51441ca0015SHong Zhang int PetscML_getrow(ML_Operator *ML_data, int N_requested_rows, int requested_rows[], 5155582bec1SHong Zhang int allocated_space, int columns[], double values[], int row_lengths[]) 5165582bec1SHong Zhang { 5175582bec1SHong Zhang PetscErrorCode ierr; 5185582bec1SHong Zhang Mat Aloc; 5195582bec1SHong Zhang Mat_SeqAIJ *a; 5205582bec1SHong Zhang PetscInt m,i,j,k=0,row,*aj; 5215582bec1SHong Zhang PetscScalar *aa; 52241ca0015SHong Zhang FineGridCtx *ml=(FineGridCtx*)ML_Get_MyGetrowData(ML_data); 5235582bec1SHong Zhang 5245582bec1SHong Zhang Aloc = ml->Aloc; 5255582bec1SHong Zhang a = (Mat_SeqAIJ*)Aloc->data; 5265582bec1SHong Zhang ierr = MatGetSize(Aloc,&m,PETSC_NULL);CHKERRQ(ierr); 5275582bec1SHong Zhang 5285582bec1SHong Zhang for (i = 0; i<N_requested_rows; i++) { 5295582bec1SHong Zhang row = requested_rows[i]; 5305582bec1SHong Zhang row_lengths[i] = a->ilen[row]; 5315582bec1SHong Zhang if (allocated_space < k+row_lengths[i]) return(0); 5325582bec1SHong Zhang if ( (row >= 0) || (row <= (m-1)) ) { 5335582bec1SHong Zhang aj = a->j + a->i[row]; 5345582bec1SHong Zhang aa = a->a + a->i[row]; 5355582bec1SHong Zhang for (j=0; j<row_lengths[i]; j++){ 5365582bec1SHong Zhang columns[k] = aj[j]; 5375582bec1SHong Zhang values[k++] = aa[j]; 5385582bec1SHong Zhang } 5395582bec1SHong Zhang } 5405582bec1SHong Zhang } 5415582bec1SHong Zhang return(1); 5425582bec1SHong Zhang } 5435582bec1SHong Zhang 54441ca0015SHong Zhang int PetscML_matvec(ML_Operator *ML_data,int in_length,double p[],int out_length,double ap[]) 5455582bec1SHong Zhang { 5465582bec1SHong Zhang PetscErrorCode ierr; 54741ca0015SHong Zhang FineGridCtx *ml=(FineGridCtx*)ML_Get_MyMatvecData(ML_data); 5485582bec1SHong Zhang Mat A=ml->A, Aloc=ml->Aloc; 5495582bec1SHong Zhang PetscMPIInt size; 5505582bec1SHong Zhang PetscScalar *pwork=ml->pwork; 5515582bec1SHong Zhang PetscInt i; 5525582bec1SHong Zhang 5537adad957SLisandro Dalcin ierr = MPI_Comm_size(((PetscObject)A)->comm,&size);CHKERRQ(ierr); 5545582bec1SHong Zhang if (size == 1){ 55524a42b14SHong Zhang ierr = VecPlaceArray(ml->x,p);CHKERRQ(ierr); 5565582bec1SHong Zhang } else { 5575582bec1SHong Zhang for (i=0; i<in_length; i++) pwork[i] = p[i]; 5585582bec1SHong Zhang PetscML_comm(pwork,ml); 55924a42b14SHong Zhang ierr = VecPlaceArray(ml->x,pwork);CHKERRQ(ierr); 5605582bec1SHong Zhang } 56124a42b14SHong Zhang ierr = VecPlaceArray(ml->y,ap);CHKERRQ(ierr); 56224a42b14SHong Zhang ierr = MatMult(Aloc,ml->x,ml->y);CHKERRQ(ierr); 563957c941cSHong Zhang ierr = VecResetArray(ml->x);CHKERRQ(ierr); 564957c941cSHong Zhang ierr = VecResetArray(ml->y);CHKERRQ(ierr); 5655582bec1SHong Zhang return 0; 5665582bec1SHong Zhang } 5675582bec1SHong Zhang 5685582bec1SHong Zhang int PetscML_comm(double p[],void *ML_data) 5695582bec1SHong Zhang { 5705582bec1SHong Zhang PetscErrorCode ierr; 5715582bec1SHong Zhang FineGridCtx *ml=(FineGridCtx*)ML_data; 5725582bec1SHong Zhang Mat A=ml->A; 5735582bec1SHong Zhang Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data; 5745582bec1SHong Zhang PetscMPIInt size; 575d0f46423SBarry Smith PetscInt i,in_length=A->rmap->n,out_length=ml->Aloc->cmap->n; 5765582bec1SHong Zhang PetscScalar *array; 5775582bec1SHong Zhang 5787adad957SLisandro Dalcin ierr = MPI_Comm_size(((PetscObject)A)->comm,&size);CHKERRQ(ierr); 5795582bec1SHong Zhang if (size == 1) return 0; 58024a42b14SHong Zhang 58124a42b14SHong Zhang ierr = VecPlaceArray(ml->y,p);CHKERRQ(ierr); 582ca9f406cSSatish Balay ierr = VecScatterBegin(a->Mvctx,ml->y,a->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 583ca9f406cSSatish Balay ierr = VecScatterEnd(a->Mvctx,ml->y,a->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 5847d15518fSHong Zhang ierr = VecResetArray(ml->y);CHKERRQ(ierr); 5855582bec1SHong Zhang ierr = VecGetArray(a->lvec,&array);CHKERRQ(ierr); 5865582bec1SHong Zhang for (i=in_length; i<out_length; i++){ 5875582bec1SHong Zhang p[i] = array[i-in_length]; 5885582bec1SHong Zhang } 5897d15518fSHong Zhang ierr = VecRestoreArray(a->lvec,&array);CHKERRQ(ierr); 5905582bec1SHong Zhang return 0; 5915582bec1SHong Zhang } 5925582bec1SHong Zhang #undef __FUNCT__ 5935582bec1SHong Zhang #define __FUNCT__ "MatMult_ML" 5945582bec1SHong Zhang PetscErrorCode MatMult_ML(Mat A,Vec x,Vec y) 5955582bec1SHong Zhang { 5965582bec1SHong Zhang PetscErrorCode ierr; 5975582bec1SHong Zhang Mat_MLShell *shell; 5985582bec1SHong Zhang PetscScalar *xarray,*yarray; 5995582bec1SHong Zhang PetscInt x_length,y_length; 6005582bec1SHong Zhang 6015582bec1SHong Zhang PetscFunctionBegin; 602a146b4dcSHong Zhang ierr = MatShellGetContext(A,(void **)&shell);CHKERRQ(ierr); 6035582bec1SHong Zhang ierr = VecGetArray(x,&xarray);CHKERRQ(ierr); 6045582bec1SHong Zhang ierr = VecGetArray(y,&yarray);CHKERRQ(ierr); 6055582bec1SHong Zhang x_length = shell->mlmat->invec_leng; 6065582bec1SHong Zhang y_length = shell->mlmat->outvec_leng; 6075582bec1SHong Zhang 6085582bec1SHong Zhang ML_Operator_Apply(shell->mlmat,x_length,xarray,y_length,yarray); 6095582bec1SHong Zhang 6105582bec1SHong Zhang ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr); 6115582bec1SHong Zhang ierr = VecRestoreArray(y,&yarray);CHKERRQ(ierr); 6125582bec1SHong Zhang PetscFunctionReturn(0); 6135582bec1SHong Zhang } 6145582bec1SHong Zhang /* MatMultAdd_ML - Compute y = w + A*x */ 6155582bec1SHong Zhang #undef __FUNCT__ 6165582bec1SHong Zhang #define __FUNCT__ "MatMultAdd_ML" 6175582bec1SHong Zhang PetscErrorCode MatMultAdd_ML(Mat A,Vec x,Vec w,Vec y) 6185582bec1SHong Zhang { 6195582bec1SHong Zhang PetscErrorCode ierr; 6205582bec1SHong Zhang Mat_MLShell *shell; 6215582bec1SHong Zhang PetscScalar *xarray,*yarray; 6225582bec1SHong Zhang PetscInt x_length,y_length; 6235582bec1SHong Zhang 6245582bec1SHong Zhang PetscFunctionBegin; 625a146b4dcSHong Zhang ierr = MatShellGetContext(A,(void **)&shell);CHKERRQ(ierr); 6265582bec1SHong Zhang ierr = VecGetArray(x,&xarray);CHKERRQ(ierr); 6275582bec1SHong Zhang ierr = VecGetArray(y,&yarray);CHKERRQ(ierr); 6285582bec1SHong Zhang 6295582bec1SHong Zhang x_length = shell->mlmat->invec_leng; 6305582bec1SHong Zhang y_length = shell->mlmat->outvec_leng; 6315582bec1SHong Zhang 6325582bec1SHong Zhang ML_Operator_Apply(shell->mlmat,x_length,xarray,y_length,yarray); 6335582bec1SHong Zhang 6345582bec1SHong Zhang ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr); 6355582bec1SHong Zhang ierr = VecRestoreArray(y,&yarray);CHKERRQ(ierr); 636efb30889SBarry Smith ierr = VecAXPY(y,1.0,w);CHKERRQ(ierr); 6375582bec1SHong Zhang 6385582bec1SHong Zhang PetscFunctionReturn(0); 6395582bec1SHong Zhang } 6405582bec1SHong Zhang 6415582bec1SHong Zhang /* newtype is ignored because "ml" is not listed under Petsc MatType yet */ 6425582bec1SHong Zhang #undef __FUNCT__ 6435582bec1SHong Zhang #define __FUNCT__ "MatConvert_MPIAIJ_ML" 64475179d2cSHong Zhang PetscErrorCode MatConvert_MPIAIJ_ML(Mat A,MatType newtype,MatReuse scall,Mat *Aloc) 6455582bec1SHong Zhang { 6465582bec1SHong Zhang PetscErrorCode ierr; 6475582bec1SHong Zhang Mat_MPIAIJ *mpimat=(Mat_MPIAIJ*)A->data; 6485582bec1SHong Zhang Mat_SeqAIJ *mat,*a=(Mat_SeqAIJ*)(mpimat->A)->data,*b=(Mat_SeqAIJ*)(mpimat->B)->data; 6495582bec1SHong Zhang PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j; 6505582bec1SHong Zhang PetscScalar *aa=a->a,*ba=b->a,*ca; 651d0f46423SBarry Smith PetscInt am=A->rmap->n,an=A->cmap->n,i,j,k; 6525582bec1SHong Zhang PetscInt *ci,*cj,ncols; 6535582bec1SHong Zhang 6545582bec1SHong Zhang PetscFunctionBegin; 6555582bec1SHong Zhang if (am != an) SETERRQ2(PETSC_ERR_ARG_WRONG,"A must have a square diagonal portion, am: %d != an: %d",am,an); 6565582bec1SHong Zhang 6575582bec1SHong Zhang if (scall == MAT_INITIAL_MATRIX){ 6585582bec1SHong Zhang ierr = PetscMalloc((1+am)*sizeof(PetscInt),&ci);CHKERRQ(ierr); 6595582bec1SHong Zhang ci[0] = 0; 6605582bec1SHong Zhang for (i=0; i<am; i++){ 6615582bec1SHong Zhang ci[i+1] = ci[i] + (ai[i+1] - ai[i]) + (bi[i+1] - bi[i]); 6625582bec1SHong Zhang } 6635582bec1SHong Zhang ierr = PetscMalloc((1+ci[am])*sizeof(PetscInt),&cj);CHKERRQ(ierr); 6645582bec1SHong Zhang ierr = PetscMalloc((1+ci[am])*sizeof(PetscScalar),&ca);CHKERRQ(ierr); 6655582bec1SHong Zhang 6665582bec1SHong Zhang k = 0; 6675582bec1SHong Zhang for (i=0; i<am; i++){ 6685582bec1SHong Zhang /* diagonal portion of A */ 6695582bec1SHong Zhang ncols = ai[i+1] - ai[i]; 6705582bec1SHong Zhang for (j=0; j<ncols; j++) { 6715582bec1SHong Zhang cj[k] = *aj++; 6725582bec1SHong Zhang ca[k++] = *aa++; 6735582bec1SHong Zhang } 6745582bec1SHong Zhang /* off-diagonal portion of A */ 6755582bec1SHong Zhang ncols = bi[i+1] - bi[i]; 6765582bec1SHong Zhang for (j=0; j<ncols; j++) { 6775582bec1SHong Zhang cj[k] = an + (*bj); bj++; 6785582bec1SHong Zhang ca[k++] = *ba++; 6795582bec1SHong Zhang } 6805582bec1SHong Zhang } 6815582bec1SHong Zhang if (k != ci[am]) SETERRQ2(PETSC_ERR_ARG_WRONG,"k: %d != ci[am]: %d",k,ci[am]); 6825582bec1SHong Zhang 6835582bec1SHong Zhang /* put together the new matrix */ 684d0f46423SBarry Smith an = mpimat->A->cmap->n+mpimat->B->cmap->n; 6855582bec1SHong Zhang ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,an,ci,cj,ca,Aloc);CHKERRQ(ierr); 6865582bec1SHong Zhang 6875582bec1SHong Zhang /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 6885582bec1SHong Zhang /* Since these are PETSc arrays, change flags to free them as necessary. */ 6895582bec1SHong Zhang mat = (Mat_SeqAIJ*)(*Aloc)->data; 6903756052fSSatish Balay mat->free_a = PETSC_TRUE; 6913756052fSSatish Balay mat->free_ij = PETSC_TRUE; 6923756052fSSatish Balay 6935582bec1SHong Zhang mat->nonew = 0; 6945582bec1SHong Zhang } else if (scall == MAT_REUSE_MATRIX){ 6955582bec1SHong Zhang mat=(Mat_SeqAIJ*)(*Aloc)->data; 6965582bec1SHong Zhang ci = mat->i; cj = mat->j; ca = mat->a; 6975582bec1SHong Zhang for (i=0; i<am; i++) { 6985582bec1SHong Zhang /* diagonal portion of A */ 6995582bec1SHong Zhang ncols = ai[i+1] - ai[i]; 7005582bec1SHong Zhang for (j=0; j<ncols; j++) *ca++ = *aa++; 7015582bec1SHong Zhang /* off-diagonal portion of A */ 7025582bec1SHong Zhang ncols = bi[i+1] - bi[i]; 7035582bec1SHong Zhang for (j=0; j<ncols; j++) *ca++ = *ba++; 7045582bec1SHong Zhang } 7055582bec1SHong Zhang } else { 7065582bec1SHong Zhang SETERRQ1(PETSC_ERR_ARG_WRONG,"Invalid MatReuse %d",(int)scall); 7075582bec1SHong Zhang } 7085582bec1SHong Zhang PetscFunctionReturn(0); 7095582bec1SHong Zhang } 7105582bec1SHong Zhang extern PetscErrorCode MatDestroy_Shell(Mat); 7115582bec1SHong Zhang #undef __FUNCT__ 7125582bec1SHong Zhang #define __FUNCT__ "MatDestroy_ML" 7135582bec1SHong Zhang PetscErrorCode MatDestroy_ML(Mat A) 7145582bec1SHong Zhang { 7155582bec1SHong Zhang PetscErrorCode ierr; 7165582bec1SHong Zhang Mat_MLShell *shell; 7175582bec1SHong Zhang 7185582bec1SHong Zhang PetscFunctionBegin; 719a146b4dcSHong Zhang ierr = MatShellGetContext(A,(void **)&shell);CHKERRQ(ierr); 7205582bec1SHong Zhang ierr = VecDestroy(shell->y);CHKERRQ(ierr); 7215582bec1SHong Zhang ierr = PetscFree(shell);CHKERRQ(ierr); 7225582bec1SHong Zhang ierr = MatDestroy_Shell(A);CHKERRQ(ierr); 723dbd8c25aSHong Zhang ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr); 7245582bec1SHong Zhang PetscFunctionReturn(0); 7255582bec1SHong Zhang } 7265582bec1SHong Zhang 7275582bec1SHong Zhang #undef __FUNCT__ 728eef31507SHong Zhang #define __FUNCT__ "MatWrapML_SeqAIJ" 729573998d7SHong Zhang PetscErrorCode MatWrapML_SeqAIJ(ML_Operator *mlmat,MatReuse reuse,Mat *newmat) 7305582bec1SHong Zhang { 731e14861a4SHong Zhang struct ML_CSR_MSRdata *matdata = (struct ML_CSR_MSRdata *)mlmat->data; 7325582bec1SHong Zhang PetscErrorCode ierr; 7333e826d7bSSatish Balay PetscInt m=mlmat->outvec_leng,n=mlmat->invec_leng,*nnz,nz_max; 734aa85bbbfSHong Zhang PetscInt *ml_cols=matdata->columns,*ml_rowptr=matdata->rowptr,*aj,i,j,k; 735e14861a4SHong Zhang PetscScalar *ml_vals=matdata->values,*aa; 7365582bec1SHong Zhang 7375582bec1SHong Zhang PetscFunctionBegin; 738e14861a4SHong Zhang if ( mlmat->getrow == NULL) SETERRQ(PETSC_ERR_ARG_NULL,"mlmat->getrow = NULL"); 739ab718edeSHong Zhang if (m != n){ /* ML Pmat and Rmat are in CSR format. Pass array pointers into SeqAIJ matrix */ 740573998d7SHong Zhang if (reuse){ 741573998d7SHong Zhang Mat_SeqAIJ *aij= (Mat_SeqAIJ*)(*newmat)->data; 742aa85bbbfSHong Zhang aij->i = ml_rowptr; 743573998d7SHong Zhang aij->j = ml_cols; 744573998d7SHong Zhang aij->a = ml_vals; 745573998d7SHong Zhang } else { 746aa85bbbfSHong Zhang /* sort ml_cols and ml_vals */ 747aa85bbbfSHong Zhang ierr = PetscMalloc((m+1)*sizeof(PetscInt),&nnz); 748aa85bbbfSHong Zhang for (i=0; i<m; i++) { 749aa85bbbfSHong Zhang nnz[i] = ml_rowptr[i+1] - ml_rowptr[i]; 750aa85bbbfSHong Zhang } 751aa85bbbfSHong Zhang aj = ml_cols; aa = ml_vals; 752aa85bbbfSHong Zhang for (i=0; i<m; i++){ 753aa85bbbfSHong Zhang ierr = PetscSortIntWithScalarArray(nnz[i],aj,aa);CHKERRQ(ierr); 754aa85bbbfSHong Zhang aj += nnz[i]; aa += nnz[i]; 755aa85bbbfSHong Zhang } 756aa85bbbfSHong Zhang ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,m,n,ml_rowptr,ml_cols,ml_vals,newmat);CHKERRQ(ierr); 757aa85bbbfSHong Zhang ierr = PetscFree(nnz);CHKERRQ(ierr); 758573998d7SHong Zhang } 75924a42b14SHong Zhang PetscFunctionReturn(0); 76024a42b14SHong Zhang } 7615582bec1SHong Zhang 762e14861a4SHong Zhang /* ML Amat is in MSR format. Copy its data into SeqAIJ matrix */ 763f69a0ea3SMatthew Knepley ierr = MatCreate(PETSC_COMM_SELF,newmat);CHKERRQ(ierr); 764f69a0ea3SMatthew Knepley ierr = MatSetSizes(*newmat,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 7655582bec1SHong Zhang ierr = MatSetType(*newmat,MATSEQAIJ);CHKERRQ(ierr); 766e14861a4SHong Zhang 767573998d7SHong Zhang ierr = PetscMalloc((m+1)*sizeof(PetscInt),&nnz); 768573998d7SHong Zhang nz_max = 1; 7695582bec1SHong Zhang for (i=0; i<m; i++) { 770e14861a4SHong Zhang nnz[i] = ml_cols[i+1] - ml_cols[i] + 1; 771573998d7SHong Zhang if (nnz[i] > nz_max) nz_max += nnz[i]; 7725582bec1SHong Zhang } 7735582bec1SHong Zhang 774573998d7SHong Zhang ierr = MatSeqAIJSetPreallocation(*newmat,0,nnz);CHKERRQ(ierr); 7757c307921SBarry Smith ierr = PetscMalloc(nz_max*(sizeof(PetscInt)+sizeof(PetscScalar)),&aj);CHKERRQ(ierr); 776e14861a4SHong Zhang aa = (PetscScalar*)(aj + nz_max); 77724a42b14SHong Zhang 7785582bec1SHong Zhang for (i=0; i<m; i++){ 779e14861a4SHong Zhang k = 0; 780e14861a4SHong Zhang /* diagonal entry */ 781e14861a4SHong Zhang aj[k] = i; aa[k++] = ml_vals[i]; 782e14861a4SHong Zhang /* off diagonal entries */ 783e14861a4SHong Zhang for (j=ml_cols[i]; j<ml_cols[i+1]; j++){ 784e14861a4SHong Zhang aj[k] = ml_cols[j]; aa[k++] = ml_vals[j]; 78524a42b14SHong Zhang } 786ab718edeSHong Zhang /* sort aj and aa */ 787e14861a4SHong Zhang ierr = PetscSortIntWithScalarArray(nnz[i],aj,aa);CHKERRQ(ierr); 788e14861a4SHong Zhang ierr = MatSetValues(*newmat,1,&i,nnz[i],aj,aa,INSERT_VALUES);CHKERRQ(ierr); 7895582bec1SHong Zhang } 7905582bec1SHong Zhang ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 7915582bec1SHong Zhang ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 792573998d7SHong Zhang 793e14861a4SHong Zhang ierr = PetscFree(aj);CHKERRQ(ierr); 7943e826d7bSSatish Balay ierr = PetscFree(nnz);CHKERRQ(ierr); 7955582bec1SHong Zhang PetscFunctionReturn(0); 7965582bec1SHong Zhang } 7975582bec1SHong Zhang 7985582bec1SHong Zhang #undef __FUNCT__ 799eef31507SHong Zhang #define __FUNCT__ "MatWrapML_SHELL" 800573998d7SHong Zhang PetscErrorCode MatWrapML_SHELL(ML_Operator *mlmat,MatReuse reuse,Mat *newmat) 8015582bec1SHong Zhang { 8025582bec1SHong Zhang PetscErrorCode ierr; 8035582bec1SHong Zhang PetscInt m,n; 8045582bec1SHong Zhang ML_Comm *MLcomm; 8055582bec1SHong Zhang Mat_MLShell *shellctx; 8065582bec1SHong Zhang 8075582bec1SHong Zhang PetscFunctionBegin; 8085582bec1SHong Zhang m = mlmat->outvec_leng; 8095582bec1SHong Zhang n = mlmat->invec_leng; 8105582bec1SHong Zhang if (!m || !n){ 8115582bec1SHong Zhang newmat = PETSC_NULL; 812573998d7SHong Zhang PetscFunctionReturn(0); 813573998d7SHong Zhang } 814573998d7SHong Zhang 815573998d7SHong Zhang if (reuse){ 816573998d7SHong Zhang ierr = MatShellGetContext(*newmat,(void **)&shellctx);CHKERRQ(ierr); 817573998d7SHong Zhang shellctx->mlmat = mlmat; 818573998d7SHong Zhang PetscFunctionReturn(0); 819573998d7SHong Zhang } 820573998d7SHong Zhang 8215582bec1SHong Zhang MLcomm = mlmat->comm; 8225582bec1SHong Zhang ierr = PetscNew(Mat_MLShell,&shellctx);CHKERRQ(ierr); 8235582bec1SHong Zhang ierr = MatCreateShell(MLcomm->USR_comm,m,n,PETSC_DETERMINE,PETSC_DETERMINE,shellctx,newmat);CHKERRQ(ierr); 8243e826d7bSSatish Balay ierr = MatShellSetOperation(*newmat,MATOP_MULT,(void(*)(void))MatMult_ML);CHKERRQ(ierr); 8253e826d7bSSatish Balay ierr = MatShellSetOperation(*newmat,MATOP_MULT_ADD,(void(*)(void))MatMultAdd_ML);CHKERRQ(ierr); 8265582bec1SHong Zhang shellctx->A = *newmat; 8275582bec1SHong Zhang shellctx->mlmat = mlmat; 8285582bec1SHong Zhang ierr = VecCreate(PETSC_COMM_WORLD,&shellctx->y);CHKERRQ(ierr); 8295582bec1SHong Zhang ierr = VecSetSizes(shellctx->y,m,PETSC_DECIDE);CHKERRQ(ierr); 8305582bec1SHong Zhang ierr = VecSetFromOptions(shellctx->y);CHKERRQ(ierr); 8315582bec1SHong Zhang (*newmat)->ops->destroy = MatDestroy_ML; 8325582bec1SHong Zhang PetscFunctionReturn(0); 8335582bec1SHong Zhang } 834e14861a4SHong Zhang 835e14861a4SHong Zhang #undef __FUNCT__ 836eef31507SHong Zhang #define __FUNCT__ "MatWrapML_MPIAIJ" 837eef31507SHong Zhang PetscErrorCode MatWrapML_MPIAIJ(ML_Operator *mlmat,Mat *newmat) 838e14861a4SHong Zhang { 839ab718edeSHong Zhang struct ML_CSR_MSRdata *matdata = (struct ML_CSR_MSRdata *)mlmat->data; 840ab718edeSHong Zhang PetscInt *ml_cols=matdata->columns,*aj; 841ab718edeSHong Zhang PetscScalar *ml_vals=matdata->values,*aa; 842e14861a4SHong Zhang PetscErrorCode ierr; 843ab718edeSHong Zhang PetscInt i,j,k,*gordering; 8443e826d7bSSatish Balay PetscInt m=mlmat->outvec_leng,n,*nnzA,*nnzB,*nnz,nz_max,row; 845ab718edeSHong Zhang Mat A; 846e14861a4SHong Zhang 847e14861a4SHong Zhang PetscFunctionBegin; 848ab718edeSHong Zhang if (mlmat->getrow == NULL) SETERRQ(PETSC_ERR_ARG_NULL,"mlmat->getrow = NULL"); 849ab718edeSHong Zhang n = mlmat->invec_leng; 850e14861a4SHong Zhang if (m != n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"m %d must equal to n %d",m,n); 851e14861a4SHong Zhang 852f69a0ea3SMatthew Knepley ierr = MatCreate(mlmat->comm->USR_comm,&A);CHKERRQ(ierr); 853f69a0ea3SMatthew Knepley ierr = MatSetSizes(A,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 854ab718edeSHong Zhang ierr = MatSetType(A,MATMPIAIJ);CHKERRQ(ierr); 8553e826d7bSSatish Balay ierr = PetscMalloc3(m,PetscInt,&nnzA,m,PetscInt,&nnzB,m,PetscInt,&nnz);CHKERRQ(ierr); 8563e826d7bSSatish Balay 857e14861a4SHong Zhang nz_max = 0; 858e14861a4SHong Zhang for (i=0; i<m; i++){ 859ab718edeSHong Zhang nnz[i] = ml_cols[i+1] - ml_cols[i] + 1; 860e14861a4SHong Zhang if (nz_max < nnz[i]) nz_max = nnz[i]; 861ab718edeSHong Zhang nnzA[i] = 1; /* diag */ 862ab718edeSHong Zhang for (j=ml_cols[i]; j<ml_cols[i+1]; j++){ 863ab718edeSHong Zhang if (ml_cols[j] < m) nnzA[i]++; 864e14861a4SHong Zhang } 865e14861a4SHong Zhang nnzB[i] = nnz[i] - nnzA[i]; 866e14861a4SHong Zhang } 867ab718edeSHong Zhang ierr = MatMPIAIJSetPreallocation(A,0,nnzA,0,nnzB);CHKERRQ(ierr); 868e14861a4SHong Zhang 869ab718edeSHong Zhang /* insert mat values -- remap row and column indices */ 870ab718edeSHong Zhang nz_max++; 8717c307921SBarry Smith ierr = PetscMalloc(nz_max*(sizeof(PetscInt)+sizeof(PetscScalar)),&aj);CHKERRQ(ierr); 872ab718edeSHong Zhang aa = (PetscScalar*)(aj + nz_max); 873510c6b62SHong Zhang /* create global row numbering for a ML_Operator */ 874510c6b62SHong Zhang ML_build_global_numbering(mlmat,&gordering,"rows"); 875e14861a4SHong Zhang for (i=0; i<m; i++){ 876e14861a4SHong Zhang row = gordering[i]; 877ab718edeSHong Zhang k = 0; 878ab718edeSHong Zhang /* diagonal entry */ 879ab718edeSHong Zhang aj[k] = row; aa[k++] = ml_vals[i]; 880ab718edeSHong Zhang /* off diagonal entries */ 881ab718edeSHong Zhang for (j=ml_cols[i]; j<ml_cols[i+1]; j++){ 882ab718edeSHong Zhang aj[k] = gordering[ml_cols[j]]; aa[k++] = ml_vals[j]; 883e14861a4SHong Zhang } 884ab718edeSHong Zhang ierr = MatSetValues(A,1,&row,nnz[i],aj,aa,INSERT_VALUES);CHKERRQ(ierr); 885e14861a4SHong Zhang } 886ab718edeSHong Zhang ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 887ab718edeSHong Zhang ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 888ab718edeSHong Zhang *newmat = A; 889e14861a4SHong Zhang 8903e826d7bSSatish Balay ierr = PetscFree3(nnzA,nnzB,nnz); 891ab718edeSHong Zhang ierr = PetscFree(aj);CHKERRQ(ierr); 892e14861a4SHong Zhang PetscFunctionReturn(0); 893e14861a4SHong Zhang } 894