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*/ 97c4f633dSBarry Smith #include "../src/ksp/pc/impls/mg/mgimpl.h" /*I "petscmg.h" I*/ 107c4f633dSBarry Smith #include "../src/mat/impls/aij/seq/aij.h" 117c4f633dSBarry 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; 103864b637dSMatthew Knepley PetscTruth isSeq, isMPI; 1045582bec1SHong Zhang 1055582bec1SHong Zhang PetscFunctionBegin; 1065582bec1SHong Zhang ierr = PetscObjectQuery((PetscObject)pc,"PC_ML",(PetscObject *)&container);CHKERRQ(ierr); 1075582bec1SHong Zhang if (container) { 108776b82aeSLisandro Dalcin ierr = PetscContainerGetPointer(container,(void **)&pc_ml);CHKERRQ(ierr); 1095582bec1SHong Zhang } else { 1105582bec1SHong Zhang SETERRQ(PETSC_ERR_ARG_NULL,"Container does not exit"); 1115582bec1SHong Zhang } 1125582bec1SHong Zhang 113573998d7SHong Zhang if (pc->setupcalled){ 114573998d7SHong Zhang if (pc->flag == SAME_NONZERO_PATTERN){ 115573998d7SHong Zhang reuse = MAT_REUSE_MATRIX; 116573998d7SHong Zhang PetscMLdata = pc_ml->PetscMLdata; 117573998d7SHong Zhang gridctx = pc_ml->gridctx; 118573998d7SHong Zhang /* ML objects cannot be reused */ 119573998d7SHong Zhang ML_Destroy(&pc_ml->ml_object); 120573998d7SHong Zhang ML_Aggregate_Destroy(&pc_ml->agg_object); 121573998d7SHong Zhang } else { 122573998d7SHong Zhang PC_ML *pc_ml_new = PETSC_NULL; 123573998d7SHong Zhang PetscContainer container_new; 12438f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_ML,&pc_ml_new);CHKERRQ(ierr); 125573998d7SHong Zhang ierr = PetscContainerCreate(PETSC_COMM_SELF,&container_new);CHKERRQ(ierr); 126573998d7SHong Zhang ierr = PetscContainerSetPointer(container_new,pc_ml_new);CHKERRQ(ierr); 127573998d7SHong Zhang ierr = PetscContainerSetUserDestroy(container_new,PetscContainerDestroy_PC_ML);CHKERRQ(ierr); 128573998d7SHong Zhang ierr = PetscObjectCompose((PetscObject)pc,"PC_ML",(PetscObject)container_new);CHKERRQ(ierr); 129573998d7SHong Zhang 130573998d7SHong Zhang ierr = PetscMemcpy(pc_ml_new,pc_ml,sizeof(PC_ML));CHKERRQ(ierr); 131573998d7SHong Zhang ierr = PetscContainerDestroy(container);CHKERRQ(ierr); 132573998d7SHong Zhang pc_ml = pc_ml_new; 133573998d7SHong Zhang } 134573998d7SHong Zhang } 135573998d7SHong Zhang 1365582bec1SHong Zhang /* setup special features of PCML */ 1375582bec1SHong Zhang /*--------------------------------*/ 1385582bec1SHong Zhang /* covert A to Aloc to be used by ML at fine grid */ 1395582bec1SHong Zhang A = pc->pmat; 1407adad957SLisandro Dalcin ierr = MPI_Comm_size(((PetscObject)A)->comm,&size);CHKERRQ(ierr); 1415582bec1SHong Zhang pc_ml->size = size; 142864b637dSMatthew Knepley ierr = PetscTypeCompare((PetscObject) A, MATSEQAIJ, &isSeq);CHKERRQ(ierr); 143864b637dSMatthew Knepley ierr = PetscTypeCompare((PetscObject) A, MATMPIAIJ, &isMPI);CHKERRQ(ierr); 144864b637dSMatthew Knepley if (isMPI){ 145573998d7SHong Zhang if (reuse) Aloc = PetscMLdata->Aloc; 146573998d7SHong Zhang ierr = MatConvert_MPIAIJ_ML(A,PETSC_NULL,reuse,&Aloc);CHKERRQ(ierr); 147864b637dSMatthew Knepley } else if (isSeq) { 1485582bec1SHong Zhang Aloc = A; 149864b637dSMatthew Knepley } else { 150864b637dSMatthew Knepley SETERRQ(PETSC_ERR_ARG_WRONG, "Invalid matrix type for ML. ML can only handle AIJ matrices."); 1515582bec1SHong Zhang } 1525582bec1SHong Zhang 1535582bec1SHong Zhang /* create and initialize struct 'PetscMLdata' */ 154573998d7SHong Zhang if (!reuse){ 15538f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,FineGridCtx,&PetscMLdata);CHKERRQ(ierr); 1565582bec1SHong Zhang pc_ml->PetscMLdata = PetscMLdata; 157d0f46423SBarry Smith ierr = PetscMalloc((Aloc->cmap->n+1)*sizeof(PetscScalar),&PetscMLdata->pwork);CHKERRQ(ierr); 1585582bec1SHong Zhang 15924a42b14SHong Zhang ierr = VecCreate(PETSC_COMM_SELF,&PetscMLdata->x);CHKERRQ(ierr); 160d0f46423SBarry Smith ierr = VecSetSizes(PetscMLdata->x,Aloc->cmap->n,Aloc->cmap->n);CHKERRQ(ierr); 16124a42b14SHong Zhang ierr = VecSetType(PetscMLdata->x,VECSEQ);CHKERRQ(ierr); 16224a42b14SHong Zhang 16324a42b14SHong Zhang ierr = VecCreate(PETSC_COMM_SELF,&PetscMLdata->y);CHKERRQ(ierr); 164d0f46423SBarry Smith ierr = VecSetSizes(PetscMLdata->y,A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr); 16524a42b14SHong Zhang ierr = VecSetType(PetscMLdata->y,VECSEQ);CHKERRQ(ierr); 166573998d7SHong Zhang } 167573998d7SHong Zhang PetscMLdata->A = A; 168573998d7SHong Zhang PetscMLdata->Aloc = Aloc; 16924a42b14SHong Zhang 1705582bec1SHong Zhang /* create ML discretization matrix at fine grid */ 17145cf47abSHong Zhang /* ML requires input of fine-grid matrix. It determines nlevels. */ 1725582bec1SHong Zhang ierr = MatGetSize(Aloc,&m,&nlocal_allcols);CHKERRQ(ierr); 1735582bec1SHong Zhang ML_Create(&ml_object,pc_ml->MaxNlevels); 174573998d7SHong Zhang pc_ml->ml_object = ml_object; 1755582bec1SHong Zhang ML_Init_Amatrix(ml_object,0,m,m,PetscMLdata); 1765582bec1SHong Zhang ML_Set_Amatrix_Getrow(ml_object,0,PetscML_getrow,PetscML_comm,nlocal_allcols); 1775582bec1SHong Zhang ML_Set_Amatrix_Matvec(ml_object,0,PetscML_matvec); 1785582bec1SHong Zhang 1795582bec1SHong Zhang /* aggregation */ 1805582bec1SHong Zhang ML_Aggregate_Create(&agg_object); 181573998d7SHong Zhang pc_ml->agg_object = agg_object; 182573998d7SHong Zhang 1835582bec1SHong Zhang ML_Aggregate_Set_MaxCoarseSize(agg_object,pc_ml->MaxCoarseSize); 1845582bec1SHong Zhang /* set options */ 1855582bec1SHong Zhang switch (pc_ml->CoarsenScheme) { 1865582bec1SHong Zhang case 1: 1875582bec1SHong Zhang ML_Aggregate_Set_CoarsenScheme_Coupled(agg_object);break; 1885582bec1SHong Zhang case 2: 1895582bec1SHong Zhang ML_Aggregate_Set_CoarsenScheme_MIS(agg_object);break; 1905582bec1SHong Zhang case 3: 1915582bec1SHong Zhang ML_Aggregate_Set_CoarsenScheme_METIS(agg_object);break; 1925582bec1SHong Zhang } 1935582bec1SHong Zhang ML_Aggregate_Set_Threshold(agg_object,pc_ml->Threshold); 1945582bec1SHong Zhang ML_Aggregate_Set_DampingFactor(agg_object,pc_ml->DampingFactor); 1955582bec1SHong Zhang if (pc_ml->SpectralNormScheme_Anorm){ 1967ffd031bSHong Zhang ML_Set_SpectralNormScheme_Anorm(ml_object); 1975582bec1SHong Zhang } 1985582bec1SHong Zhang 1995582bec1SHong Zhang Nlevels = ML_Gen_MGHierarchy_UsingAggregation(ml_object,0,ML_INCREASING,agg_object); 2005582bec1SHong Zhang if (Nlevels<=0) SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Nlevels %d must > 0",Nlevels); 201573998d7SHong 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); 202573998d7SHong Zhang pc_ml->Nlevels = Nlevels; 203aa85bbbfSHong Zhang fine_level = Nlevels - 1; 204573998d7SHong Zhang if (!pc->setupcalled){ 20597177400SBarry Smith ierr = PCMGSetLevels(pc,Nlevels,PETSC_NULL);CHKERRQ(ierr); 206aa85bbbfSHong Zhang /* set default smoothers */ 207aa85bbbfSHong Zhang KSP smoother; 208aa85bbbfSHong Zhang PC subpc; 209aa85bbbfSHong Zhang for (level=1; level<=fine_level; level++){ 210aa85bbbfSHong Zhang if (size == 1){ 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 } else { 216aa85bbbfSHong Zhang ierr = PCMGGetSmoother(pc,level,&smoother);CHKERRQ(ierr); 217aa85bbbfSHong Zhang ierr = KSPSetType(smoother,KSPRICHARDSON);CHKERRQ(ierr); 218aa85bbbfSHong Zhang ierr = KSPGetPC(smoother,&subpc);CHKERRQ(ierr); 219aa85bbbfSHong Zhang ierr = PCSetType(subpc,PCSOR);CHKERRQ(ierr); 220aa85bbbfSHong Zhang } 221aa85bbbfSHong Zhang } 22297177400SBarry Smith ierr = PCSetFromOptions_MG(pc);CHKERRQ(ierr); /* should be called in PCSetFromOptions_ML(), but cannot be called prior to PCMGSetLevels() */ 223573998d7SHong Zhang } 2245582bec1SHong Zhang 225573998d7SHong Zhang if (!reuse){ 2265582bec1SHong Zhang ierr = PetscMalloc(Nlevels*sizeof(GridCtx),&gridctx);CHKERRQ(ierr); 2275582bec1SHong Zhang pc_ml->gridctx = gridctx; 228573998d7SHong Zhang } 2295582bec1SHong Zhang 2305582bec1SHong Zhang /* wrap ML matrices by PETSc shell matrices at coarsened grids. 2315582bec1SHong Zhang Level 0 is the finest grid for ML, but coarsest for PETSc! */ 232e14861a4SHong Zhang gridctx[fine_level].A = A; 233573998d7SHong Zhang 234e14861a4SHong Zhang level = fine_level - 1; 235ab718edeSHong Zhang if (size == 1){ /* convert ML P, R and A into seqaij format */ 2365582bec1SHong Zhang for (mllevel=1; mllevel<Nlevels; mllevel++){ 237e14861a4SHong Zhang mlmat = &(ml_object->Pmat[mllevel]); 238573998d7SHong Zhang ierr = MatWrapML_SeqAIJ(mlmat,reuse,&gridctx[level].P);CHKERRQ(ierr); 239e14861a4SHong Zhang mlmat = &(ml_object->Rmat[mllevel-1]); 240573998d7SHong Zhang ierr = MatWrapML_SeqAIJ(mlmat,reuse,&gridctx[level].R);CHKERRQ(ierr); 241573998d7SHong Zhang 242573998d7SHong Zhang mlmat = &(ml_object->Amat[mllevel]); 243573998d7SHong Zhang if (reuse){ 244573998d7SHong Zhang /* ML matrix A changes sparse pattern although PETSc A doesn't, thus gridctx[level].A must be recreated! */ 245573998d7SHong Zhang ierr = MatDestroy(gridctx[level].A);CHKERRQ(ierr); 246573998d7SHong Zhang } 247573998d7SHong Zhang ierr = MatWrapML_SeqAIJ(mlmat,MAT_INITIAL_MATRIX,&gridctx[level].A);CHKERRQ(ierr); 2485582bec1SHong Zhang level--; 2495582bec1SHong Zhang } 250ab718edeSHong Zhang } else { /* convert ML P and R into shell format, ML A into mpiaij format */ 2515582bec1SHong Zhang for (mllevel=1; mllevel<Nlevels; mllevel++){ 2525582bec1SHong Zhang mlmat = &(ml_object->Pmat[mllevel]); 253573998d7SHong Zhang ierr = MatWrapML_SHELL(mlmat,reuse,&gridctx[level].P);CHKERRQ(ierr); 254ab718edeSHong Zhang mlmat = &(ml_object->Rmat[mllevel-1]); 255573998d7SHong Zhang ierr = MatWrapML_SHELL(mlmat,reuse,&gridctx[level].R);CHKERRQ(ierr); 256573998d7SHong Zhang 2575582bec1SHong Zhang mlmat = &(ml_object->Amat[mllevel]); 258573998d7SHong Zhang if (reuse){ 259573998d7SHong Zhang ierr = MatDestroy(gridctx[level].A);CHKERRQ(ierr); 260573998d7SHong Zhang } 261eef31507SHong Zhang ierr = MatWrapML_MPIAIJ(mlmat,&gridctx[level].A);CHKERRQ(ierr); 2625582bec1SHong Zhang level--; 2635582bec1SHong Zhang } 2645582bec1SHong Zhang } 2655582bec1SHong Zhang 266573998d7SHong Zhang /* create vectors and ksp at all levels */ 267573998d7SHong Zhang if (!reuse){ 268ac346b81SHong Zhang for (level=0; level<fine_level; level++){ 269573998d7SHong Zhang level1 = level + 1; 270e64afeacSLisandro Dalcin ierr = VecCreate(((PetscObject)gridctx[level].A)->comm,&gridctx[level].x);CHKERRQ(ierr); 271d0f46423SBarry Smith ierr = VecSetSizes(gridctx[level].x,gridctx[level].A->cmap->n,PETSC_DECIDE);CHKERRQ(ierr); 2725582bec1SHong Zhang ierr = VecSetType(gridctx[level].x,VECMPI);CHKERRQ(ierr); 27397177400SBarry Smith ierr = PCMGSetX(pc,level,gridctx[level].x);CHKERRQ(ierr); 2745582bec1SHong Zhang 275e64afeacSLisandro Dalcin ierr = VecCreate(((PetscObject)gridctx[level].A)->comm,&gridctx[level].b);CHKERRQ(ierr); 276d0f46423SBarry Smith ierr = VecSetSizes(gridctx[level].b,gridctx[level].A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr); 2775582bec1SHong Zhang ierr = VecSetType(gridctx[level].b,VECMPI);CHKERRQ(ierr); 27897177400SBarry Smith ierr = PCMGSetRhs(pc,level,gridctx[level].b);CHKERRQ(ierr); 279ac346b81SHong Zhang 280e64afeacSLisandro Dalcin ierr = VecCreate(((PetscObject)gridctx[level1].A)->comm,&gridctx[level1].r);CHKERRQ(ierr); 281d0f46423SBarry Smith ierr = VecSetSizes(gridctx[level1].r,gridctx[level1].A->rmap->n,PETSC_DECIDE);CHKERRQ(ierr); 282ac346b81SHong Zhang ierr = VecSetType(gridctx[level1].r,VECMPI);CHKERRQ(ierr); 28397177400SBarry Smith ierr = PCMGSetR(pc,level1,gridctx[level1].r);CHKERRQ(ierr); 284ac346b81SHong Zhang 2855582bec1SHong Zhang if (level == 0){ 28697177400SBarry Smith ierr = PCMGGetCoarseSolve(pc,&gridctx[level].ksp);CHKERRQ(ierr); 2875582bec1SHong Zhang } else { 28897177400SBarry Smith ierr = PCMGGetSmoother(pc,level,&gridctx[level].ksp);CHKERRQ(ierr); 289573998d7SHong Zhang } 290573998d7SHong Zhang } 291573998d7SHong Zhang ierr = PCMGGetSmoother(pc,fine_level,&gridctx[fine_level].ksp);CHKERRQ(ierr); 292573998d7SHong Zhang } 293573998d7SHong Zhang 294573998d7SHong Zhang /* create coarse level and the interpolation between the levels */ 295573998d7SHong Zhang for (level=0; level<fine_level; level++){ 296573998d7SHong Zhang level1 = level + 1; 297aea2a34eSBarry Smith ierr = PCMGSetInterpolation(pc,level1,gridctx[level].P);CHKERRQ(ierr); 298573998d7SHong Zhang ierr = PCMGSetRestriction(pc,level1,gridctx[level].R);CHKERRQ(ierr); 299573998d7SHong Zhang if (level > 0){ 30097177400SBarry Smith ierr = PCMGSetResidual(pc,level,PCMGDefaultResidual,gridctx[level].A);CHKERRQ(ierr); 3015582bec1SHong Zhang } 3025582bec1SHong Zhang ierr = KSPSetOperators(gridctx[level].ksp,gridctx[level].A,gridctx[level].A,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 3035582bec1SHong Zhang } 30497177400SBarry Smith ierr = PCMGSetResidual(pc,fine_level,PCMGDefaultResidual,gridctx[fine_level].A);CHKERRQ(ierr); 305ac346b81SHong Zhang ierr = KSPSetOperators(gridctx[fine_level].ksp,gridctx[level].A,gridctx[fine_level].A,DIFFERENT_NONZERO_PATTERN);CHKERRQ(ierr); 3065582bec1SHong Zhang 3075582bec1SHong Zhang /* now call PCSetUp_MG() */ 308573998d7SHong Zhang /*-------------------------------*/ 3095582bec1SHong Zhang ierr = (*pc_ml->PCSetUp)(pc);CHKERRQ(ierr); 3105582bec1SHong Zhang PetscFunctionReturn(0); 3115582bec1SHong Zhang } 3125582bec1SHong Zhang 3135582bec1SHong Zhang #undef __FUNCT__ 314776b82aeSLisandro Dalcin #define __FUNCT__ "PetscContainerDestroy_PC_ML" 315776b82aeSLisandro Dalcin PetscErrorCode PetscContainerDestroy_PC_ML(void *ptr) 3165582bec1SHong Zhang { 3175582bec1SHong Zhang PetscErrorCode ierr; 3185582bec1SHong Zhang PC_ML *pc_ml = (PC_ML*)ptr; 319573998d7SHong Zhang PetscInt level,fine_level=pc_ml->Nlevels-1; 3205582bec1SHong Zhang 3215582bec1SHong Zhang PetscFunctionBegin; 3225582bec1SHong Zhang ML_Aggregate_Destroy(&pc_ml->agg_object); 3235582bec1SHong Zhang ML_Destroy(&pc_ml->ml_object); 3245582bec1SHong Zhang 32538f2d2fdSLisandro Dalcin if (pc_ml->PetscMLdata) { 3265582bec1SHong Zhang ierr = PetscFree(pc_ml->PetscMLdata->pwork);CHKERRQ(ierr); 32738f2d2fdSLisandro Dalcin if (pc_ml->size > 1) {ierr = MatDestroy(pc_ml->PetscMLdata->Aloc);CHKERRQ(ierr);} 328ac346b81SHong Zhang if (pc_ml->PetscMLdata->x){ierr = VecDestroy(pc_ml->PetscMLdata->x);CHKERRQ(ierr);} 329ac346b81SHong Zhang if (pc_ml->PetscMLdata->y){ierr = VecDestroy(pc_ml->PetscMLdata->y);CHKERRQ(ierr);} 33038f2d2fdSLisandro Dalcin } 3315582bec1SHong Zhang ierr = PetscFree(pc_ml->PetscMLdata);CHKERRQ(ierr); 3325582bec1SHong Zhang 333573998d7SHong Zhang for (level=0; level<fine_level; level++){ 334ac346b81SHong Zhang if (pc_ml->gridctx[level].A){ierr = MatDestroy(pc_ml->gridctx[level].A);CHKERRQ(ierr);} 335ac346b81SHong Zhang if (pc_ml->gridctx[level].P){ierr = MatDestroy(pc_ml->gridctx[level].P);CHKERRQ(ierr);} 336ac346b81SHong Zhang if (pc_ml->gridctx[level].R){ierr = MatDestroy(pc_ml->gridctx[level].R);CHKERRQ(ierr);} 337ac346b81SHong Zhang if (pc_ml->gridctx[level].x){ierr = VecDestroy(pc_ml->gridctx[level].x);CHKERRQ(ierr);} 338ac346b81SHong Zhang if (pc_ml->gridctx[level].b){ierr = VecDestroy(pc_ml->gridctx[level].b);CHKERRQ(ierr);} 339ac346b81SHong Zhang if (pc_ml->gridctx[level+1].r){ierr = VecDestroy(pc_ml->gridctx[level+1].r);CHKERRQ(ierr);} 3405582bec1SHong Zhang } 3415582bec1SHong Zhang ierr = PetscFree(pc_ml->gridctx);CHKERRQ(ierr); 3425582bec1SHong Zhang ierr = PetscFree(pc_ml);CHKERRQ(ierr); 3435582bec1SHong Zhang PetscFunctionReturn(0); 3445582bec1SHong Zhang } 3455582bec1SHong Zhang /* -------------------------------------------------------------------------- */ 3465582bec1SHong Zhang /* 3475582bec1SHong Zhang PCDestroy_ML - Destroys the private context for the ML preconditioner 3485582bec1SHong Zhang that was created with PCCreate_ML(). 3495582bec1SHong Zhang 3505582bec1SHong Zhang Input Parameter: 3515582bec1SHong Zhang . pc - the preconditioner context 3525582bec1SHong Zhang 3535582bec1SHong Zhang Application Interface Routine: PCDestroy() 3545582bec1SHong Zhang */ 3555582bec1SHong Zhang #undef __FUNCT__ 3565582bec1SHong Zhang #define __FUNCT__ "PCDestroy_ML" 3576ca4d86aSHong Zhang PetscErrorCode PCDestroy_ML(PC pc) 3585582bec1SHong Zhang { 3595582bec1SHong Zhang PetscErrorCode ierr; 3605582bec1SHong Zhang PC_ML *pc_ml=PETSC_NULL; 361776b82aeSLisandro Dalcin PetscContainer container; 3625582bec1SHong Zhang 3635582bec1SHong Zhang PetscFunctionBegin; 3645582bec1SHong Zhang ierr = PetscObjectQuery((PetscObject)pc,"PC_ML",(PetscObject *)&container);CHKERRQ(ierr); 3655582bec1SHong Zhang if (container) { 366776b82aeSLisandro Dalcin ierr = PetscContainerGetPointer(container,(void **)&pc_ml);CHKERRQ(ierr); 3675582bec1SHong Zhang pc->ops->destroy = pc_ml->PCDestroy; 3685582bec1SHong Zhang } else { 3695582bec1SHong Zhang SETERRQ(PETSC_ERR_ARG_NULL,"Container does not exit"); 3705582bec1SHong Zhang } 3719cb0ec93SHong Zhang /* detach pc and PC_ML and dereference container */ 37238f2d2fdSLisandro Dalcin ierr = PetscContainerDestroy(container);CHKERRQ(ierr); 3735582bec1SHong Zhang ierr = PetscObjectCompose((PetscObject)pc,"PC_ML",0);CHKERRQ(ierr); 37438f2d2fdSLisandro Dalcin if (pc->ops->destroy) { 3755582bec1SHong Zhang ierr = (*pc->ops->destroy)(pc);CHKERRQ(ierr); 37638f2d2fdSLisandro Dalcin } 3775582bec1SHong Zhang PetscFunctionReturn(0); 3785582bec1SHong Zhang } 3795582bec1SHong Zhang 3805582bec1SHong Zhang #undef __FUNCT__ 3815582bec1SHong Zhang #define __FUNCT__ "PCSetFromOptions_ML" 3826ca4d86aSHong Zhang PetscErrorCode PCSetFromOptions_ML(PC pc) 3835582bec1SHong Zhang { 3845582bec1SHong Zhang PetscErrorCode ierr; 38538f2d2fdSLisandro Dalcin PetscInt indx,m,PrintLevel; 3865582bec1SHong Zhang PetscTruth flg; 3875582bec1SHong Zhang const char *scheme[] = {"Uncoupled","Coupled","MIS","METIS"}; 3885582bec1SHong Zhang PC_ML *pc_ml=PETSC_NULL; 389776b82aeSLisandro Dalcin PetscContainer container; 3909dcbbd2bSBarry Smith PCMGType mgtype; 3915582bec1SHong Zhang 3925582bec1SHong Zhang PetscFunctionBegin; 3935582bec1SHong Zhang ierr = PetscObjectQuery((PetscObject)pc,"PC_ML",(PetscObject *)&container);CHKERRQ(ierr); 3945582bec1SHong Zhang if (container) { 395776b82aeSLisandro Dalcin ierr = PetscContainerGetPointer(container,(void **)&pc_ml);CHKERRQ(ierr); 3965582bec1SHong Zhang } else { 3975582bec1SHong Zhang SETERRQ(PETSC_ERR_ARG_NULL,"Container does not exit"); 3985582bec1SHong Zhang } 3996ca4d86aSHong Zhang 4005582bec1SHong Zhang /* inherited MG options */ 4016ca4d86aSHong Zhang ierr = PetscOptionsHead("Multigrid options(inherited)");CHKERRQ(ierr); 4025582bec1SHong Zhang ierr = PetscOptionsInt("-pc_mg_cycles","1 for V cycle, 2 for W-cycle","MGSetCycles",1,&m,&flg);CHKERRQ(ierr); 4035582bec1SHong Zhang ierr = PetscOptionsInt("-pc_mg_smoothup","Number of post-smoothing steps","MGSetNumberSmoothUp",1,&m,&flg);CHKERRQ(ierr); 4045582bec1SHong Zhang ierr = PetscOptionsInt("-pc_mg_smoothdown","Number of pre-smoothing steps","MGSetNumberSmoothDown",1,&m,&flg);CHKERRQ(ierr); 4059dcbbd2bSBarry Smith ierr = PetscOptionsEnum("-pc_mg_type","Multigrid type","PCMGSetType",PCMGTypes,(PetscEnum)PC_MG_MULTIPLICATIVE,(PetscEnum*)&mgtype,&flg);CHKERRQ(ierr); 4065582bec1SHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 4075582bec1SHong Zhang 4085582bec1SHong Zhang /* ML options */ 4095582bec1SHong Zhang ierr = PetscOptionsHead("ML options");CHKERRQ(ierr); 4105582bec1SHong Zhang /* set defaults */ 4115582bec1SHong Zhang PrintLevel = 0; 4125582bec1SHong Zhang indx = 0; 4135582bec1SHong Zhang ierr = PetscOptionsInt("-pc_ml_PrintLevel","Print level","ML_Set_PrintLevel",PrintLevel,&PrintLevel,PETSC_NULL);CHKERRQ(ierr); 4145582bec1SHong Zhang ML_Set_PrintLevel(PrintLevel); 415573998d7SHong Zhang ierr = PetscOptionsInt("-pc_ml_maxNlevels","Maximum number of levels","None",pc_ml->MaxNlevels,&pc_ml->MaxNlevels,PETSC_NULL);CHKERRQ(ierr); 416573998d7SHong Zhang ierr = PetscOptionsInt("-pc_ml_maxCoarseSize","Maximum coarsest mesh size","ML_Aggregate_Set_MaxCoarseSize",pc_ml->MaxCoarseSize,&pc_ml->MaxCoarseSize,PETSC_NULL);CHKERRQ(ierr); 417573998d7SHong Zhang ierr = PetscOptionsEList("-pc_ml_CoarsenScheme","Aggregate Coarsen Scheme","ML_Aggregate_Set_CoarsenScheme_*",scheme,4,scheme[0],&indx,PETSC_NULL);CHKERRQ(ierr); /* ??? */ 4185582bec1SHong Zhang pc_ml->CoarsenScheme = indx; 4195582bec1SHong Zhang 420573998d7SHong Zhang ierr = PetscOptionsReal("-pc_ml_DampingFactor","P damping factor","ML_Aggregate_Set_DampingFactor",pc_ml->DampingFactor,&pc_ml->DampingFactor,PETSC_NULL);CHKERRQ(ierr); 4215582bec1SHong Zhang 422573998d7SHong Zhang ierr = PetscOptionsReal("-pc_ml_Threshold","Smoother drop tol","ML_Aggregate_Set_Threshold",pc_ml->Threshold,&pc_ml->Threshold,PETSC_NULL);CHKERRQ(ierr); 4235582bec1SHong Zhang 42440597110SHong 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); 4255582bec1SHong Zhang 4265582bec1SHong Zhang ierr = PetscOptionsTail();CHKERRQ(ierr); 4275582bec1SHong Zhang PetscFunctionReturn(0); 4285582bec1SHong Zhang } 4295582bec1SHong Zhang 4305582bec1SHong Zhang /* -------------------------------------------------------------------------- */ 4315582bec1SHong Zhang /* 4325582bec1SHong Zhang PCCreate_ML - Creates a ML preconditioner context, PC_ML, 4335582bec1SHong Zhang and sets this as the private data within the generic preconditioning 4345582bec1SHong Zhang context, PC, that was created within PCCreate(). 4355582bec1SHong Zhang 4365582bec1SHong Zhang Input Parameter: 4375582bec1SHong Zhang . pc - the preconditioner context 4385582bec1SHong Zhang 4395582bec1SHong Zhang Application Interface Routine: PCCreate() 4405582bec1SHong Zhang */ 4415582bec1SHong Zhang 4425582bec1SHong Zhang /*MC 4431e5ab15bSHong Zhang PCML - Use algebraic multigrid preconditioning. This preconditioner requires you provide 4445582bec1SHong Zhang fine grid discretization matrix. The coarser grid matrices and restriction/interpolation 4456ca4d86aSHong Zhang operators are computed by ML, with the matrices coverted to PETSc matrices in aij format 4466ca4d86aSHong Zhang and the restriction/interpolation operators wrapped as PETSc shell matrices. 4475582bec1SHong Zhang 4486ca4d86aSHong Zhang Options Database Key: 4496ca4d86aSHong Zhang Multigrid options(inherited) 4506ca4d86aSHong Zhang + -pc_mg_cycles <1>: 1 for V cycle, 2 for W-cycle (MGSetCycles) 4516ca4d86aSHong Zhang . -pc_mg_smoothup <1>: Number of post-smoothing steps (MGSetNumberSmoothUp) 4526ca4d86aSHong Zhang . -pc_mg_smoothdown <1>: Number of pre-smoothing steps (MGSetNumberSmoothDown) 453f41ab451SVictor Eijkhout - -pc_mg_type <multiplicative>: (one of) additive multiplicative full cascade kascade 4546ca4d86aSHong Zhang 45551f519a2SBarry Smith ML options: 4566ca4d86aSHong Zhang + -pc_ml_PrintLevel <0>: Print level (ML_Set_PrintLevel) 4576ca4d86aSHong Zhang . -pc_ml_maxNlevels <10>: Maximum number of levels (None) 4586ca4d86aSHong Zhang . -pc_ml_maxCoarseSize <1>: Maximum coarsest mesh size (ML_Aggregate_Set_MaxCoarseSize) 459f41ab451SVictor Eijkhout . -pc_ml_CoarsenScheme <Uncoupled>: (one of) Uncoupled Coupled MIS METIS 4606ca4d86aSHong Zhang . -pc_ml_DampingFactor <1.33333>: P damping factor (ML_Aggregate_Set_DampingFactor) 4616ca4d86aSHong Zhang . -pc_ml_Threshold <0>: Smoother drop tol (ML_Aggregate_Set_Threshold) 4627ffd031bSHong Zhang - -pc_ml_SpectralNormScheme_Anorm <false>: Method used for estimating spectral radius (ML_Set_SpectralNormScheme_Anorm) 4635582bec1SHong Zhang 4645582bec1SHong Zhang Level: intermediate 4655582bec1SHong Zhang 4665582bec1SHong Zhang Concepts: multigrid 4675582bec1SHong Zhang 4685582bec1SHong Zhang .seealso: PCCreate(), PCSetType(), PCType (for list of available types), PC, PCMGType, 46997177400SBarry Smith PCMGSetLevels(), PCMGGetLevels(), PCMGSetType(), MPSetCycles(), PCMGSetNumberSmoothDown(), 47097177400SBarry Smith PCMGSetNumberSmoothUp(), PCMGGetCoarseSolve(), PCMGSetResidual(), PCMGSetInterpolation(), 47197177400SBarry Smith PCMGSetRestriction(), PCMGGetSmoother(), PCMGGetSmootherUp(), PCMGGetSmootherDown(), 47297177400SBarry Smith PCMGSetCyclesOnLevel(), PCMGSetRhs(), PCMGSetX(), PCMGSetR() 4735582bec1SHong Zhang M*/ 4745582bec1SHong Zhang 4755582bec1SHong Zhang EXTERN_C_BEGIN 4765582bec1SHong Zhang #undef __FUNCT__ 4775582bec1SHong Zhang #define __FUNCT__ "PCCreate_ML" 478dba47a55SKris Buschelman PetscErrorCode PETSCKSP_DLLEXPORT PCCreate_ML(PC pc) 4795582bec1SHong Zhang { 4805582bec1SHong Zhang PetscErrorCode ierr; 4815582bec1SHong Zhang PC_ML *pc_ml; 482776b82aeSLisandro Dalcin PetscContainer container; 4835582bec1SHong Zhang 4845582bec1SHong Zhang PetscFunctionBegin; 485573998d7SHong Zhang /* PCML is an inherited class of PCMG. Initialize pc as PCMG */ 4865582bec1SHong Zhang ierr = PCSetType(pc,PCMG);CHKERRQ(ierr); /* calls PCCreate_MG() and MGCreate_Private() */ 4875582bec1SHong Zhang 4885582bec1SHong Zhang /* create a supporting struct and attach it to pc */ 48938f2d2fdSLisandro Dalcin ierr = PetscNewLog(pc,PC_ML,&pc_ml);CHKERRQ(ierr); 490776b82aeSLisandro Dalcin ierr = PetscContainerCreate(PETSC_COMM_SELF,&container);CHKERRQ(ierr); 491776b82aeSLisandro Dalcin ierr = PetscContainerSetPointer(container,pc_ml);CHKERRQ(ierr); 492776b82aeSLisandro Dalcin ierr = PetscContainerSetUserDestroy(container,PetscContainerDestroy_PC_ML);CHKERRQ(ierr); 4935582bec1SHong Zhang ierr = PetscObjectCompose((PetscObject)pc,"PC_ML",(PetscObject)container);CHKERRQ(ierr); 4945582bec1SHong Zhang 495573998d7SHong Zhang pc_ml->ml_object = 0; 496573998d7SHong Zhang pc_ml->agg_object = 0; 497573998d7SHong Zhang pc_ml->gridctx = 0; 498573998d7SHong Zhang pc_ml->PetscMLdata = 0; 499573998d7SHong Zhang pc_ml->Nlevels = -1; 500573998d7SHong Zhang pc_ml->MaxNlevels = 10; 501573998d7SHong Zhang pc_ml->MaxCoarseSize = 1; 502573998d7SHong Zhang pc_ml->CoarsenScheme = 1; /* ??? */ 503573998d7SHong Zhang pc_ml->Threshold = 0.0; 504573998d7SHong Zhang pc_ml->DampingFactor = 4.0/3.0; 505573998d7SHong Zhang pc_ml->SpectralNormScheme_Anorm = PETSC_FALSE; 506573998d7SHong Zhang pc_ml->size = 0; 507573998d7SHong Zhang 5085582bec1SHong Zhang pc_ml->PCSetUp = pc->ops->setup; 5095582bec1SHong Zhang pc_ml->PCDestroy = pc->ops->destroy; 5105582bec1SHong Zhang 5115582bec1SHong Zhang /* overwrite the pointers of PCMG by the functions of PCML */ 5125582bec1SHong Zhang pc->ops->setfromoptions = PCSetFromOptions_ML; 5135582bec1SHong Zhang pc->ops->setup = PCSetUp_ML; 5145582bec1SHong Zhang pc->ops->destroy = PCDestroy_ML; 5155582bec1SHong Zhang PetscFunctionReturn(0); 5165582bec1SHong Zhang } 5175582bec1SHong Zhang EXTERN_C_END 5185582bec1SHong Zhang 51941ca0015SHong Zhang int PetscML_getrow(ML_Operator *ML_data, int N_requested_rows, int requested_rows[], 5205582bec1SHong Zhang int allocated_space, int columns[], double values[], int row_lengths[]) 5215582bec1SHong Zhang { 5225582bec1SHong Zhang PetscErrorCode ierr; 5235582bec1SHong Zhang Mat Aloc; 5245582bec1SHong Zhang Mat_SeqAIJ *a; 5255582bec1SHong Zhang PetscInt m,i,j,k=0,row,*aj; 5265582bec1SHong Zhang PetscScalar *aa; 52741ca0015SHong Zhang FineGridCtx *ml=(FineGridCtx*)ML_Get_MyGetrowData(ML_data); 5285582bec1SHong Zhang 5295582bec1SHong Zhang Aloc = ml->Aloc; 5305582bec1SHong Zhang a = (Mat_SeqAIJ*)Aloc->data; 5315582bec1SHong Zhang ierr = MatGetSize(Aloc,&m,PETSC_NULL);CHKERRQ(ierr); 5325582bec1SHong Zhang 5335582bec1SHong Zhang for (i = 0; i<N_requested_rows; i++) { 5345582bec1SHong Zhang row = requested_rows[i]; 5355582bec1SHong Zhang row_lengths[i] = a->ilen[row]; 5365582bec1SHong Zhang if (allocated_space < k+row_lengths[i]) return(0); 5375582bec1SHong Zhang if ( (row >= 0) || (row <= (m-1)) ) { 5385582bec1SHong Zhang aj = a->j + a->i[row]; 5395582bec1SHong Zhang aa = a->a + a->i[row]; 5405582bec1SHong Zhang for (j=0; j<row_lengths[i]; j++){ 5415582bec1SHong Zhang columns[k] = aj[j]; 5425582bec1SHong Zhang values[k++] = aa[j]; 5435582bec1SHong Zhang } 5445582bec1SHong Zhang } 5455582bec1SHong Zhang } 5465582bec1SHong Zhang return(1); 5475582bec1SHong Zhang } 5485582bec1SHong Zhang 54941ca0015SHong Zhang int PetscML_matvec(ML_Operator *ML_data,int in_length,double p[],int out_length,double ap[]) 5505582bec1SHong Zhang { 5515582bec1SHong Zhang PetscErrorCode ierr; 55241ca0015SHong Zhang FineGridCtx *ml=(FineGridCtx*)ML_Get_MyMatvecData(ML_data); 5535582bec1SHong Zhang Mat A=ml->A, Aloc=ml->Aloc; 5545582bec1SHong Zhang PetscMPIInt size; 5555582bec1SHong Zhang PetscScalar *pwork=ml->pwork; 5565582bec1SHong Zhang PetscInt i; 5575582bec1SHong Zhang 5587adad957SLisandro Dalcin ierr = MPI_Comm_size(((PetscObject)A)->comm,&size);CHKERRQ(ierr); 5595582bec1SHong Zhang if (size == 1){ 56024a42b14SHong Zhang ierr = VecPlaceArray(ml->x,p);CHKERRQ(ierr); 5615582bec1SHong Zhang } else { 5625582bec1SHong Zhang for (i=0; i<in_length; i++) pwork[i] = p[i]; 5635582bec1SHong Zhang PetscML_comm(pwork,ml); 56424a42b14SHong Zhang ierr = VecPlaceArray(ml->x,pwork);CHKERRQ(ierr); 5655582bec1SHong Zhang } 56624a42b14SHong Zhang ierr = VecPlaceArray(ml->y,ap);CHKERRQ(ierr); 56724a42b14SHong Zhang ierr = MatMult(Aloc,ml->x,ml->y);CHKERRQ(ierr); 568957c941cSHong Zhang ierr = VecResetArray(ml->x);CHKERRQ(ierr); 569957c941cSHong Zhang ierr = VecResetArray(ml->y);CHKERRQ(ierr); 5705582bec1SHong Zhang return 0; 5715582bec1SHong Zhang } 5725582bec1SHong Zhang 5735582bec1SHong Zhang int PetscML_comm(double p[],void *ML_data) 5745582bec1SHong Zhang { 5755582bec1SHong Zhang PetscErrorCode ierr; 5765582bec1SHong Zhang FineGridCtx *ml=(FineGridCtx*)ML_data; 5775582bec1SHong Zhang Mat A=ml->A; 5785582bec1SHong Zhang Mat_MPIAIJ *a = (Mat_MPIAIJ*)A->data; 5795582bec1SHong Zhang PetscMPIInt size; 580d0f46423SBarry Smith PetscInt i,in_length=A->rmap->n,out_length=ml->Aloc->cmap->n; 5815582bec1SHong Zhang PetscScalar *array; 5825582bec1SHong Zhang 5837adad957SLisandro Dalcin ierr = MPI_Comm_size(((PetscObject)A)->comm,&size);CHKERRQ(ierr); 5845582bec1SHong Zhang if (size == 1) return 0; 58524a42b14SHong Zhang 58624a42b14SHong Zhang ierr = VecPlaceArray(ml->y,p);CHKERRQ(ierr); 587ca9f406cSSatish Balay ierr = VecScatterBegin(a->Mvctx,ml->y,a->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 588ca9f406cSSatish Balay ierr = VecScatterEnd(a->Mvctx,ml->y,a->lvec,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 5897d15518fSHong Zhang ierr = VecResetArray(ml->y);CHKERRQ(ierr); 5905582bec1SHong Zhang ierr = VecGetArray(a->lvec,&array);CHKERRQ(ierr); 5915582bec1SHong Zhang for (i=in_length; i<out_length; i++){ 5925582bec1SHong Zhang p[i] = array[i-in_length]; 5935582bec1SHong Zhang } 5947d15518fSHong Zhang ierr = VecRestoreArray(a->lvec,&array);CHKERRQ(ierr); 5955582bec1SHong Zhang return 0; 5965582bec1SHong Zhang } 5975582bec1SHong Zhang #undef __FUNCT__ 5985582bec1SHong Zhang #define __FUNCT__ "MatMult_ML" 5995582bec1SHong Zhang PetscErrorCode MatMult_ML(Mat A,Vec x,Vec y) 6005582bec1SHong Zhang { 6015582bec1SHong Zhang PetscErrorCode ierr; 6025582bec1SHong Zhang Mat_MLShell *shell; 6035582bec1SHong Zhang PetscScalar *xarray,*yarray; 6045582bec1SHong Zhang PetscInt x_length,y_length; 6055582bec1SHong Zhang 6065582bec1SHong Zhang PetscFunctionBegin; 607a146b4dcSHong Zhang ierr = MatShellGetContext(A,(void **)&shell);CHKERRQ(ierr); 6085582bec1SHong Zhang ierr = VecGetArray(x,&xarray);CHKERRQ(ierr); 6095582bec1SHong Zhang ierr = VecGetArray(y,&yarray);CHKERRQ(ierr); 6105582bec1SHong Zhang x_length = shell->mlmat->invec_leng; 6115582bec1SHong Zhang y_length = shell->mlmat->outvec_leng; 6125582bec1SHong Zhang 6135582bec1SHong Zhang ML_Operator_Apply(shell->mlmat,x_length,xarray,y_length,yarray); 6145582bec1SHong Zhang 6155582bec1SHong Zhang ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr); 6165582bec1SHong Zhang ierr = VecRestoreArray(y,&yarray);CHKERRQ(ierr); 6175582bec1SHong Zhang PetscFunctionReturn(0); 6185582bec1SHong Zhang } 6195582bec1SHong Zhang /* MatMultAdd_ML - Compute y = w + A*x */ 6205582bec1SHong Zhang #undef __FUNCT__ 6215582bec1SHong Zhang #define __FUNCT__ "MatMultAdd_ML" 6225582bec1SHong Zhang PetscErrorCode MatMultAdd_ML(Mat A,Vec x,Vec w,Vec y) 6235582bec1SHong Zhang { 6245582bec1SHong Zhang PetscErrorCode ierr; 6255582bec1SHong Zhang Mat_MLShell *shell; 6265582bec1SHong Zhang PetscScalar *xarray,*yarray; 6275582bec1SHong Zhang PetscInt x_length,y_length; 6285582bec1SHong Zhang 6295582bec1SHong Zhang PetscFunctionBegin; 630a146b4dcSHong Zhang ierr = MatShellGetContext(A,(void **)&shell);CHKERRQ(ierr); 6315582bec1SHong Zhang ierr = VecGetArray(x,&xarray);CHKERRQ(ierr); 6325582bec1SHong Zhang ierr = VecGetArray(y,&yarray);CHKERRQ(ierr); 6335582bec1SHong Zhang 6345582bec1SHong Zhang x_length = shell->mlmat->invec_leng; 6355582bec1SHong Zhang y_length = shell->mlmat->outvec_leng; 6365582bec1SHong Zhang 6375582bec1SHong Zhang ML_Operator_Apply(shell->mlmat,x_length,xarray,y_length,yarray); 6385582bec1SHong Zhang 6395582bec1SHong Zhang ierr = VecRestoreArray(x,&xarray);CHKERRQ(ierr); 6405582bec1SHong Zhang ierr = VecRestoreArray(y,&yarray);CHKERRQ(ierr); 641efb30889SBarry Smith ierr = VecAXPY(y,1.0,w);CHKERRQ(ierr); 6425582bec1SHong Zhang 6435582bec1SHong Zhang PetscFunctionReturn(0); 6445582bec1SHong Zhang } 6455582bec1SHong Zhang 6465582bec1SHong Zhang /* newtype is ignored because "ml" is not listed under Petsc MatType yet */ 6475582bec1SHong Zhang #undef __FUNCT__ 6485582bec1SHong Zhang #define __FUNCT__ "MatConvert_MPIAIJ_ML" 64975179d2cSHong Zhang PetscErrorCode MatConvert_MPIAIJ_ML(Mat A,MatType newtype,MatReuse scall,Mat *Aloc) 6505582bec1SHong Zhang { 6515582bec1SHong Zhang PetscErrorCode ierr; 6525582bec1SHong Zhang Mat_MPIAIJ *mpimat=(Mat_MPIAIJ*)A->data; 6535582bec1SHong Zhang Mat_SeqAIJ *mat,*a=(Mat_SeqAIJ*)(mpimat->A)->data,*b=(Mat_SeqAIJ*)(mpimat->B)->data; 6545582bec1SHong Zhang PetscInt *ai=a->i,*aj=a->j,*bi=b->i,*bj=b->j; 6555582bec1SHong Zhang PetscScalar *aa=a->a,*ba=b->a,*ca; 656d0f46423SBarry Smith PetscInt am=A->rmap->n,an=A->cmap->n,i,j,k; 6575582bec1SHong Zhang PetscInt *ci,*cj,ncols; 6585582bec1SHong Zhang 6595582bec1SHong Zhang PetscFunctionBegin; 6605582bec1SHong Zhang if (am != an) SETERRQ2(PETSC_ERR_ARG_WRONG,"A must have a square diagonal portion, am: %d != an: %d",am,an); 6615582bec1SHong Zhang 6625582bec1SHong Zhang if (scall == MAT_INITIAL_MATRIX){ 6635582bec1SHong Zhang ierr = PetscMalloc((1+am)*sizeof(PetscInt),&ci);CHKERRQ(ierr); 6645582bec1SHong Zhang ci[0] = 0; 6655582bec1SHong Zhang for (i=0; i<am; i++){ 6665582bec1SHong Zhang ci[i+1] = ci[i] + (ai[i+1] - ai[i]) + (bi[i+1] - bi[i]); 6675582bec1SHong Zhang } 6685582bec1SHong Zhang ierr = PetscMalloc((1+ci[am])*sizeof(PetscInt),&cj);CHKERRQ(ierr); 6695582bec1SHong Zhang ierr = PetscMalloc((1+ci[am])*sizeof(PetscScalar),&ca);CHKERRQ(ierr); 6705582bec1SHong Zhang 6715582bec1SHong Zhang k = 0; 6725582bec1SHong Zhang for (i=0; i<am; i++){ 6735582bec1SHong Zhang /* diagonal portion of A */ 6745582bec1SHong Zhang ncols = ai[i+1] - ai[i]; 6755582bec1SHong Zhang for (j=0; j<ncols; j++) { 6765582bec1SHong Zhang cj[k] = *aj++; 6775582bec1SHong Zhang ca[k++] = *aa++; 6785582bec1SHong Zhang } 6795582bec1SHong Zhang /* off-diagonal portion of A */ 6805582bec1SHong Zhang ncols = bi[i+1] - bi[i]; 6815582bec1SHong Zhang for (j=0; j<ncols; j++) { 6825582bec1SHong Zhang cj[k] = an + (*bj); bj++; 6835582bec1SHong Zhang ca[k++] = *ba++; 6845582bec1SHong Zhang } 6855582bec1SHong Zhang } 6865582bec1SHong Zhang if (k != ci[am]) SETERRQ2(PETSC_ERR_ARG_WRONG,"k: %d != ci[am]: %d",k,ci[am]); 6875582bec1SHong Zhang 6885582bec1SHong Zhang /* put together the new matrix */ 689d0f46423SBarry Smith an = mpimat->A->cmap->n+mpimat->B->cmap->n; 6905582bec1SHong Zhang ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,am,an,ci,cj,ca,Aloc);CHKERRQ(ierr); 6915582bec1SHong Zhang 6925582bec1SHong Zhang /* MatCreateSeqAIJWithArrays flags matrix so PETSc doesn't free the user's arrays. */ 6935582bec1SHong Zhang /* Since these are PETSc arrays, change flags to free them as necessary. */ 6945582bec1SHong Zhang mat = (Mat_SeqAIJ*)(*Aloc)->data; 6953756052fSSatish Balay mat->free_a = PETSC_TRUE; 6963756052fSSatish Balay mat->free_ij = PETSC_TRUE; 6973756052fSSatish Balay 6985582bec1SHong Zhang mat->nonew = 0; 6995582bec1SHong Zhang } else if (scall == MAT_REUSE_MATRIX){ 7005582bec1SHong Zhang mat=(Mat_SeqAIJ*)(*Aloc)->data; 7015582bec1SHong Zhang ci = mat->i; cj = mat->j; ca = mat->a; 7025582bec1SHong Zhang for (i=0; i<am; i++) { 7035582bec1SHong Zhang /* diagonal portion of A */ 7045582bec1SHong Zhang ncols = ai[i+1] - ai[i]; 7055582bec1SHong Zhang for (j=0; j<ncols; j++) *ca++ = *aa++; 7065582bec1SHong Zhang /* off-diagonal portion of A */ 7075582bec1SHong Zhang ncols = bi[i+1] - bi[i]; 7085582bec1SHong Zhang for (j=0; j<ncols; j++) *ca++ = *ba++; 7095582bec1SHong Zhang } 7105582bec1SHong Zhang } else { 7115582bec1SHong Zhang SETERRQ1(PETSC_ERR_ARG_WRONG,"Invalid MatReuse %d",(int)scall); 7125582bec1SHong Zhang } 7135582bec1SHong Zhang PetscFunctionReturn(0); 7145582bec1SHong Zhang } 7155582bec1SHong Zhang extern PetscErrorCode MatDestroy_Shell(Mat); 7165582bec1SHong Zhang #undef __FUNCT__ 7175582bec1SHong Zhang #define __FUNCT__ "MatDestroy_ML" 7185582bec1SHong Zhang PetscErrorCode MatDestroy_ML(Mat A) 7195582bec1SHong Zhang { 7205582bec1SHong Zhang PetscErrorCode ierr; 7215582bec1SHong Zhang Mat_MLShell *shell; 7225582bec1SHong Zhang 7235582bec1SHong Zhang PetscFunctionBegin; 724a146b4dcSHong Zhang ierr = MatShellGetContext(A,(void **)&shell);CHKERRQ(ierr); 7255582bec1SHong Zhang ierr = VecDestroy(shell->y);CHKERRQ(ierr); 7265582bec1SHong Zhang ierr = PetscFree(shell);CHKERRQ(ierr); 7275582bec1SHong Zhang ierr = MatDestroy_Shell(A);CHKERRQ(ierr); 728dbd8c25aSHong Zhang ierr = PetscObjectChangeTypeName((PetscObject)A,0);CHKERRQ(ierr); 7295582bec1SHong Zhang PetscFunctionReturn(0); 7305582bec1SHong Zhang } 7315582bec1SHong Zhang 7325582bec1SHong Zhang #undef __FUNCT__ 733eef31507SHong Zhang #define __FUNCT__ "MatWrapML_SeqAIJ" 734573998d7SHong Zhang PetscErrorCode MatWrapML_SeqAIJ(ML_Operator *mlmat,MatReuse reuse,Mat *newmat) 7355582bec1SHong Zhang { 736e14861a4SHong Zhang struct ML_CSR_MSRdata *matdata = (struct ML_CSR_MSRdata *)mlmat->data; 7375582bec1SHong Zhang PetscErrorCode ierr; 7383e826d7bSSatish Balay PetscInt m=mlmat->outvec_leng,n=mlmat->invec_leng,*nnz,nz_max; 739aa85bbbfSHong Zhang PetscInt *ml_cols=matdata->columns,*ml_rowptr=matdata->rowptr,*aj,i,j,k; 740e14861a4SHong Zhang PetscScalar *ml_vals=matdata->values,*aa; 7415582bec1SHong Zhang 7425582bec1SHong Zhang PetscFunctionBegin; 743e14861a4SHong Zhang if ( mlmat->getrow == NULL) SETERRQ(PETSC_ERR_ARG_NULL,"mlmat->getrow = NULL"); 744ab718edeSHong Zhang if (m != n){ /* ML Pmat and Rmat are in CSR format. Pass array pointers into SeqAIJ matrix */ 745573998d7SHong Zhang if (reuse){ 746573998d7SHong Zhang Mat_SeqAIJ *aij= (Mat_SeqAIJ*)(*newmat)->data; 747aa85bbbfSHong Zhang aij->i = ml_rowptr; 748573998d7SHong Zhang aij->j = ml_cols; 749573998d7SHong Zhang aij->a = ml_vals; 750573998d7SHong Zhang } else { 751aa85bbbfSHong Zhang /* sort ml_cols and ml_vals */ 752aa85bbbfSHong Zhang ierr = PetscMalloc((m+1)*sizeof(PetscInt),&nnz); 753aa85bbbfSHong Zhang for (i=0; i<m; i++) { 754aa85bbbfSHong Zhang nnz[i] = ml_rowptr[i+1] - ml_rowptr[i]; 755aa85bbbfSHong Zhang } 756aa85bbbfSHong Zhang aj = ml_cols; aa = ml_vals; 757aa85bbbfSHong Zhang for (i=0; i<m; i++){ 758aa85bbbfSHong Zhang ierr = PetscSortIntWithScalarArray(nnz[i],aj,aa);CHKERRQ(ierr); 759aa85bbbfSHong Zhang aj += nnz[i]; aa += nnz[i]; 760aa85bbbfSHong Zhang } 761aa85bbbfSHong Zhang ierr = MatCreateSeqAIJWithArrays(PETSC_COMM_SELF,m,n,ml_rowptr,ml_cols,ml_vals,newmat);CHKERRQ(ierr); 762aa85bbbfSHong Zhang ierr = PetscFree(nnz);CHKERRQ(ierr); 763573998d7SHong Zhang } 76424a42b14SHong Zhang PetscFunctionReturn(0); 76524a42b14SHong Zhang } 7665582bec1SHong Zhang 767e14861a4SHong Zhang /* ML Amat is in MSR format. Copy its data into SeqAIJ matrix */ 768f69a0ea3SMatthew Knepley ierr = MatCreate(PETSC_COMM_SELF,newmat);CHKERRQ(ierr); 769f69a0ea3SMatthew Knepley ierr = MatSetSizes(*newmat,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 7705582bec1SHong Zhang ierr = MatSetType(*newmat,MATSEQAIJ);CHKERRQ(ierr); 771e14861a4SHong Zhang 772573998d7SHong Zhang ierr = PetscMalloc((m+1)*sizeof(PetscInt),&nnz); 773573998d7SHong Zhang nz_max = 1; 7745582bec1SHong Zhang for (i=0; i<m; i++) { 775e14861a4SHong Zhang nnz[i] = ml_cols[i+1] - ml_cols[i] + 1; 776573998d7SHong Zhang if (nnz[i] > nz_max) nz_max += nnz[i]; 7775582bec1SHong Zhang } 7785582bec1SHong Zhang 779573998d7SHong Zhang ierr = MatSeqAIJSetPreallocation(*newmat,0,nnz);CHKERRQ(ierr); 780*1d79065fSBarry Smith ierr = PetscMalloc2(nz_max,PetscScalar,&aa,nz_max,PetscInt,&aj);CHKERRQ(ierr); 7815582bec1SHong Zhang for (i=0; i<m; i++){ 782e14861a4SHong Zhang k = 0; 783e14861a4SHong Zhang /* diagonal entry */ 784e14861a4SHong Zhang aj[k] = i; aa[k++] = ml_vals[i]; 785e14861a4SHong Zhang /* off diagonal entries */ 786e14861a4SHong Zhang for (j=ml_cols[i]; j<ml_cols[i+1]; j++){ 787e14861a4SHong Zhang aj[k] = ml_cols[j]; aa[k++] = ml_vals[j]; 78824a42b14SHong Zhang } 789ab718edeSHong Zhang /* sort aj and aa */ 790e14861a4SHong Zhang ierr = PetscSortIntWithScalarArray(nnz[i],aj,aa);CHKERRQ(ierr); 791e14861a4SHong Zhang ierr = MatSetValues(*newmat,1,&i,nnz[i],aj,aa,INSERT_VALUES);CHKERRQ(ierr); 7925582bec1SHong Zhang } 7935582bec1SHong Zhang ierr = MatAssemblyBegin(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 7945582bec1SHong Zhang ierr = MatAssemblyEnd(*newmat,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 795573998d7SHong Zhang 796*1d79065fSBarry Smith ierr = PetscFree2(aa,aj);CHKERRQ(ierr); 7973e826d7bSSatish Balay ierr = PetscFree(nnz);CHKERRQ(ierr); 7985582bec1SHong Zhang PetscFunctionReturn(0); 7995582bec1SHong Zhang } 8005582bec1SHong Zhang 8015582bec1SHong Zhang #undef __FUNCT__ 802eef31507SHong Zhang #define __FUNCT__ "MatWrapML_SHELL" 803573998d7SHong Zhang PetscErrorCode MatWrapML_SHELL(ML_Operator *mlmat,MatReuse reuse,Mat *newmat) 8045582bec1SHong Zhang { 8055582bec1SHong Zhang PetscErrorCode ierr; 8065582bec1SHong Zhang PetscInt m,n; 8075582bec1SHong Zhang ML_Comm *MLcomm; 8085582bec1SHong Zhang Mat_MLShell *shellctx; 8095582bec1SHong Zhang 8105582bec1SHong Zhang PetscFunctionBegin; 8115582bec1SHong Zhang m = mlmat->outvec_leng; 8125582bec1SHong Zhang n = mlmat->invec_leng; 8135582bec1SHong Zhang if (!m || !n){ 8145582bec1SHong Zhang newmat = PETSC_NULL; 815573998d7SHong Zhang PetscFunctionReturn(0); 816573998d7SHong Zhang } 817573998d7SHong Zhang 818573998d7SHong Zhang if (reuse){ 819573998d7SHong Zhang ierr = MatShellGetContext(*newmat,(void **)&shellctx);CHKERRQ(ierr); 820573998d7SHong Zhang shellctx->mlmat = mlmat; 821573998d7SHong Zhang PetscFunctionReturn(0); 822573998d7SHong Zhang } 823573998d7SHong Zhang 8245582bec1SHong Zhang MLcomm = mlmat->comm; 8255582bec1SHong Zhang ierr = PetscNew(Mat_MLShell,&shellctx);CHKERRQ(ierr); 8265582bec1SHong Zhang ierr = MatCreateShell(MLcomm->USR_comm,m,n,PETSC_DETERMINE,PETSC_DETERMINE,shellctx,newmat);CHKERRQ(ierr); 8273e826d7bSSatish Balay ierr = MatShellSetOperation(*newmat,MATOP_MULT,(void(*)(void))MatMult_ML);CHKERRQ(ierr); 8283e826d7bSSatish Balay ierr = MatShellSetOperation(*newmat,MATOP_MULT_ADD,(void(*)(void))MatMultAdd_ML);CHKERRQ(ierr); 8295582bec1SHong Zhang shellctx->A = *newmat; 8305582bec1SHong Zhang shellctx->mlmat = mlmat; 8315582bec1SHong Zhang ierr = VecCreate(PETSC_COMM_WORLD,&shellctx->y);CHKERRQ(ierr); 8325582bec1SHong Zhang ierr = VecSetSizes(shellctx->y,m,PETSC_DECIDE);CHKERRQ(ierr); 8335582bec1SHong Zhang ierr = VecSetFromOptions(shellctx->y);CHKERRQ(ierr); 8345582bec1SHong Zhang (*newmat)->ops->destroy = MatDestroy_ML; 8355582bec1SHong Zhang PetscFunctionReturn(0); 8365582bec1SHong Zhang } 837e14861a4SHong Zhang 838e14861a4SHong Zhang #undef __FUNCT__ 839eef31507SHong Zhang #define __FUNCT__ "MatWrapML_MPIAIJ" 840eef31507SHong Zhang PetscErrorCode MatWrapML_MPIAIJ(ML_Operator *mlmat,Mat *newmat) 841e14861a4SHong Zhang { 842ab718edeSHong Zhang struct ML_CSR_MSRdata *matdata = (struct ML_CSR_MSRdata *)mlmat->data; 843ab718edeSHong Zhang PetscInt *ml_cols=matdata->columns,*aj; 844ab718edeSHong Zhang PetscScalar *ml_vals=matdata->values,*aa; 845e14861a4SHong Zhang PetscErrorCode ierr; 846ab718edeSHong Zhang PetscInt i,j,k,*gordering; 8473e826d7bSSatish Balay PetscInt m=mlmat->outvec_leng,n,*nnzA,*nnzB,*nnz,nz_max,row; 848ab718edeSHong Zhang Mat A; 849e14861a4SHong Zhang 850e14861a4SHong Zhang PetscFunctionBegin; 851ab718edeSHong Zhang if (mlmat->getrow == NULL) SETERRQ(PETSC_ERR_ARG_NULL,"mlmat->getrow = NULL"); 852ab718edeSHong Zhang n = mlmat->invec_leng; 853e14861a4SHong Zhang if (m != n) SETERRQ2(PETSC_ERR_ARG_OUTOFRANGE,"m %d must equal to n %d",m,n); 854e14861a4SHong Zhang 855f69a0ea3SMatthew Knepley ierr = MatCreate(mlmat->comm->USR_comm,&A);CHKERRQ(ierr); 856f69a0ea3SMatthew Knepley ierr = MatSetSizes(A,m,n,PETSC_DECIDE,PETSC_DECIDE);CHKERRQ(ierr); 857ab718edeSHong Zhang ierr = MatSetType(A,MATMPIAIJ);CHKERRQ(ierr); 8583e826d7bSSatish Balay ierr = PetscMalloc3(m,PetscInt,&nnzA,m,PetscInt,&nnzB,m,PetscInt,&nnz);CHKERRQ(ierr); 8593e826d7bSSatish Balay 860e14861a4SHong Zhang nz_max = 0; 861e14861a4SHong Zhang for (i=0; i<m; i++){ 862ab718edeSHong Zhang nnz[i] = ml_cols[i+1] - ml_cols[i] + 1; 863e14861a4SHong Zhang if (nz_max < nnz[i]) nz_max = nnz[i]; 864ab718edeSHong Zhang nnzA[i] = 1; /* diag */ 865ab718edeSHong Zhang for (j=ml_cols[i]; j<ml_cols[i+1]; j++){ 866ab718edeSHong Zhang if (ml_cols[j] < m) nnzA[i]++; 867e14861a4SHong Zhang } 868e14861a4SHong Zhang nnzB[i] = nnz[i] - nnzA[i]; 869e14861a4SHong Zhang } 870ab718edeSHong Zhang ierr = MatMPIAIJSetPreallocation(A,0,nnzA,0,nnzB);CHKERRQ(ierr); 871e14861a4SHong Zhang 872ab718edeSHong Zhang /* insert mat values -- remap row and column indices */ 873ab718edeSHong Zhang nz_max++; 874*1d79065fSBarry Smith ierr = PetscMalloc2(nz_max,PetscScalar,&aa,nz_max,PetscInt,&aj);CHKERRQ(ierr); 875510c6b62SHong Zhang /* create global row numbering for a ML_Operator */ 876510c6b62SHong Zhang ML_build_global_numbering(mlmat,&gordering,"rows"); 877e14861a4SHong Zhang for (i=0; i<m; i++){ 878e14861a4SHong Zhang row = gordering[i]; 879ab718edeSHong Zhang k = 0; 880ab718edeSHong Zhang /* diagonal entry */ 881ab718edeSHong Zhang aj[k] = row; aa[k++] = ml_vals[i]; 882ab718edeSHong Zhang /* off diagonal entries */ 883ab718edeSHong Zhang for (j=ml_cols[i]; j<ml_cols[i+1]; j++){ 884ab718edeSHong Zhang aj[k] = gordering[ml_cols[j]]; aa[k++] = ml_vals[j]; 885e14861a4SHong Zhang } 886ab718edeSHong Zhang ierr = MatSetValues(A,1,&row,nnz[i],aj,aa,INSERT_VALUES);CHKERRQ(ierr); 887e14861a4SHong Zhang } 888ab718edeSHong Zhang ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 889ab718edeSHong Zhang ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 890ab718edeSHong Zhang *newmat = A; 891e14861a4SHong Zhang 8923e826d7bSSatish Balay ierr = PetscFree3(nnzA,nnzB,nnz); 893*1d79065fSBarry Smith ierr = PetscFree2(aa,aj);CHKERRQ(ierr); 894e14861a4SHong Zhang PetscFunctionReturn(0); 895e14861a4SHong Zhang } 896